xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision dadcf80911fb48939c55327431ae8d7e47dbe367)
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;
145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm, &section));
155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalSection(dm, &sectionGlobal));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject) section), &adjSec));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(section,&pStart,&pEnd));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetChart(adjSec,pStart,pEnd));
1976185916SToby Isaac 
205f80ce2aSJacob Faibussowitsch   CHKERRQ(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 */
285f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionCreate(PetscObjectComm((PetscObject)aSec),&inverseSec));
295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetChart(inverseSec,pStart,pEnd));
305f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetLocalSize(aIS, &aSize));
315f80ce2aSJacob Faibussowitsch     CHKERRQ(ISGetIndices(aIS, &anchors));
3276185916SToby Isaac 
3376185916SToby Isaac     for (p = 0; p < aSize; p++) {
3476185916SToby Isaac       PetscInt a = anchors[p];
3576185916SToby Isaac 
365f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionAddDof(inverseSec,a,1));
3776185916SToby Isaac     }
385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(inverseSec));
395f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetStorageSize(inverseSec,&iSize));
405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(iSize,&inverse));
415f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscCalloc1(pEnd-pStart,&offsets));
425f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetChart(aSec,&aStart,&aEnd));
4376185916SToby Isaac     for (p = aStart; p < aEnd; p++) {
4476185916SToby Isaac       PetscInt dof, off;
4576185916SToby Isaac 
465f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(aSec, p, &dof));
475f80ce2aSJacob Faibussowitsch       CHKERRQ(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];
535f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(inverseSec, a, &iOff));
5476185916SToby Isaac         inverse[iOff + offsets[a-pStart]++] = p;
5576185916SToby Isaac       }
5676185916SToby Isaac     }
575f80ce2aSJacob Faibussowitsch     CHKERRQ(ISRestoreIndices(aIS, &anchors));
585f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
805f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(inverseSec,p,&iDof));
8176185916SToby Isaac       if (!iDof) continue;
825f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(inverseSec,p,&iOff));
835f80ce2aSJacob Faibussowitsch       CHKERRQ(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];
885f80ce2aSJacob Faibussowitsch         CHKERRQ(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;
965f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSectionGetDof(section,qAdj,&qAdjDof));
975f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof));
9876185916SToby Isaac           iNew += qAdjDof - qAdjCDof;
9976185916SToby Isaac         }
1005f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionAddDof(adjSec,p,iNew));
10176185916SToby Isaac       }
10276185916SToby Isaac     }
10376185916SToby Isaac 
1045f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionSetUp(adjSec));
1055f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetStorageSize(adjSec,&adjSize));
1065f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
1115f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(inverseSec,p,&iDof));
11276185916SToby Isaac       if (!iDof) continue;
1135f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(inverseSec,p,&iOff));
1145f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP));
1155f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(adjSec,p,&aDof));
1165f80ce2aSJacob Faibussowitsch       CHKERRQ(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];
1225f80ce2aSJacob Faibussowitsch         CHKERRQ(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;
1305f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSectionGetDof(section,qAdj,&qAdjDof));
1315f80ce2aSJacob Faibussowitsch           CHKERRQ(PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof));
1325f80ce2aSJacob Faibussowitsch           CHKERRQ(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       }
1385f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]));
1395f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(adjSec,p,aDof));
14076185916SToby Isaac     }
14176185916SToby Isaac     *anchorAdj = adj;
14276185916SToby Isaac 
14376185916SToby Isaac     /* clean up */
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionDestroy(&inverseSec));
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(inverse));
1465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(tmpAdjP));
1475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(tmpAdjQ));
14876185916SToby Isaac   }
14976185916SToby Isaac   else {
15076185916SToby Isaac     *anchorAdj = NULL;
1515f80ce2aSJacob Faibussowitsch     CHKERRQ(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;
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL));
1745f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetPointSF(dm, &sf));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm, &section));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetGlobalSection(dm, &sectionGlobal));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
18047a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
1815f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPIU_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm));
182cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(section, &pStart, &pEnd));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetStorageSize(section, &numDof));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &leafSectionAdj));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetChart(leafSectionAdj, 0, numDof));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &rootSectionAdj));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetChart(rootSectionAdj, 0, numDof));
189cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(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   */
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
2245f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
2255f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(section, p, &off));
2265f80ce2aSJacob Faibussowitsch     CHKERRQ(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;
2325f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(section, padj, &ndof));
2335f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetConstraintDof(section, padj, &ncdof));
234cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
2355f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof));
236cb1e1211SMatthew G Knepley       }
237cb1e1211SMatthew G Knepley     }
2385f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
23976185916SToby Isaac     if (anDof) {
24076185916SToby Isaac       for (d = off; d < off+dof; ++d) {
2415f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionAddDof(leafSectionAdj, d, anDof));
24276185916SToby Isaac       }
24376185916SToby Isaac     }
244cb1e1211SMatthew G Knepley   }
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(leafSectionAdj));
246cb1e1211SMatthew G Knepley   if (debug) {
2475f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n"));
2485f80ce2aSJacob Faibussowitsch     CHKERRQ(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) {
2525f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
2535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
254cb1e1211SMatthew G Knepley   }
255cb1e1211SMatthew G Knepley   if (debug) {
2565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n"));
2575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionView(rootSectionAdj, NULL));
258cb1e1211SMatthew G Knepley   }
259cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
260cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
26176185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
262cb1e1211SMatthew G Knepley 
2635f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
2645f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(section, p, &off));
265cb1e1211SMatthew G Knepley     if (!dof) continue;
2665f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(rootSectionAdj, off, &adof));
267cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
2685f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
269cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
270cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
271cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
272cb1e1211SMatthew G Knepley 
273cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
2745f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(section, padj, &ndof));
2755f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetConstraintDof(section, padj, &ncdof));
276cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
2775f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof));
278cb1e1211SMatthew G Knepley       }
279cb1e1211SMatthew G Knepley     }
2805f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
28176185916SToby Isaac     if (anDof) {
28276185916SToby Isaac       for (d = off; d < off+dof; ++d) {
2835f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionAddDof(rootSectionAdj, d, anDof));
28476185916SToby Isaac       }
28576185916SToby Isaac     }
286cb1e1211SMatthew G Knepley   }
2875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(rootSectionAdj));
288cb1e1211SMatthew G Knepley   if (debug) {
2895f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n"));
2905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionView(rootSectionAdj, NULL));
291cb1e1211SMatthew G Knepley   }
292cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
2935f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets));
2945f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj));
2955f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(remoteOffsets));
2966cded2eaSMatthew G. Knepley   if (debug && size > 1) {
2975f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Adjacency SF for Preallocation:\n"));
2985f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFView(sfAdj, NULL));
299cb1e1211SMatthew G Knepley   }
300cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
3015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(leafSectionAdj));
3025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetStorageSize(leafSectionAdj, &adjSize));
3035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(adjSize, &adj));
304cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
30576185916SToby Isaac     PetscInt dof, off, d, q, anDof, anOff;
30670034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
307cb1e1211SMatthew G Knepley 
308fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
3095f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
3105f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(section, p, &off));
3115f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
314cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
315cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
316cb1e1211SMatthew G Knepley 
3175f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(leafSectionAdj, d, &aoff));
318cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
319cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
320cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
321cb1e1211SMatthew G Knepley 
322cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
3235f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetDof(section, padj, &ndof));
3245f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetConstraintDof(section, padj, &ncdof));
3255f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
326cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
327cb1e1211SMatthew G Knepley           adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd;
328cb1e1211SMatthew G Knepley           ++i;
329cb1e1211SMatthew G Knepley         }
330cb1e1211SMatthew G Knepley       }
33176185916SToby Isaac       for (q = 0; q < anDof; q++) {
33276185916SToby Isaac         adj[aoff+i] = anchorAdj[anOff+q];
33376185916SToby Isaac         ++i;
33476185916SToby Isaac       }
335cb1e1211SMatthew G Knepley     }
336cb1e1211SMatthew G Knepley   }
337cb1e1211SMatthew G Knepley   /* Debugging */
338cb1e1211SMatthew G Knepley   if (debug) {
339cb1e1211SMatthew G Knepley     IS tmp;
3405f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Leaf adjacency indices\n"));
3415f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp));
3425f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(tmp, NULL));
3435f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&tmp));
344cb1e1211SMatthew G Knepley   }
345543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
3465f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetStorageSize(rootSectionAdj, &adjSize));
3475f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(adjSize, &rootAdj));
348cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
34947a6104aSMatthew G. Knepley   if (doComm) {
350543482b8SMatthew G. Knepley     const PetscInt *indegree;
351543482b8SMatthew G. Knepley     PetscInt       *remoteadj, radjsize = 0;
352543482b8SMatthew G. Knepley 
3535f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFComputeDegreeBegin(sfAdj, &indegree));
3545f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFComputeDegreeEnd(sfAdj, &indegree));
355543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
3565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscMalloc1(radjsize, &remoteadj));
3575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj));
3585f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj));
3596dba2905SMatthew G. Knepley     for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p-1])) {
360543482b8SMatthew G. Knepley       PetscInt s;
3616dba2905SMatthew G. Knepley       for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l+s] = remoteadj[r];
362543482b8SMatthew G. Knepley     }
3632c71b3e2SJacob Faibussowitsch     PetscCheckFalse(r != radjsize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", r, radjsize);
3642c71b3e2SJacob Faibussowitsch     PetscCheckFalse(l != adjSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", l, adjSize);
3655f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscFree(remoteadj));
366cb1e1211SMatthew G Knepley   }
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfAdj));
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(adj));
369cb1e1211SMatthew G Knepley   /* Debugging */
370cb1e1211SMatthew G Knepley   if (debug) {
371cb1e1211SMatthew G Knepley     IS tmp;
3725f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Root adjacency indices after gather\n"));
3735f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
3745f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(tmp, NULL));
3755f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&tmp));
376cb1e1211SMatthew G Knepley   }
377cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
378cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
37976185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
380cb1e1211SMatthew G Knepley 
3815f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
3825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(section, p, &off));
383cb1e1211SMatthew G Knepley     if (!dof) continue;
3845f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(rootSectionAdj, off, &adof));
385cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
3865f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
3875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
3885f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
389cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
390cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
391cb1e1211SMatthew G Knepley 
3925f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(rootSectionAdj, d, &adof));
3935f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
394cb1e1211SMatthew G Knepley       i    = adof-1;
39576185916SToby Isaac       for (q = 0; q < anDof; q++) {
39676185916SToby Isaac         rootAdj[aoff+i] = anchorAdj[anOff+q];
39776185916SToby Isaac         --i;
39876185916SToby Isaac       }
399cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
400cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
401cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
402cb1e1211SMatthew G Knepley 
403cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
4045f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetDof(section, padj, &ndof));
4055f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetConstraintDof(section, padj, &ncdof));
4065f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
407cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
408cb1e1211SMatthew G Knepley           rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
409cb1e1211SMatthew G Knepley           --i;
410cb1e1211SMatthew G Knepley         }
411cb1e1211SMatthew G Knepley       }
412cb1e1211SMatthew G Knepley     }
413cb1e1211SMatthew G Knepley   }
414cb1e1211SMatthew G Knepley   /* Debugging */
415cb1e1211SMatthew G Knepley   if (debug) {
416cb1e1211SMatthew G Knepley     IS tmp;
4175f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Root adjacency indices\n"));
4185f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4195f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(tmp, NULL));
4205f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&tmp));
421cb1e1211SMatthew G Knepley   }
422cb1e1211SMatthew G Knepley   /* Compress indices */
4235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(rootSectionAdj));
424cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
425cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
426cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
427cb1e1211SMatthew G Knepley 
4285f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
4295f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetConstraintDof(section, p, &cdof));
4305f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(section, p, &off));
431cb1e1211SMatthew G Knepley     if (!dof) continue;
4325f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(rootSectionAdj, off, &adof));
433cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
434cb1e1211SMatthew G Knepley     for (d = off; d < off+dof-cdof; ++d) {
4355f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(rootSectionAdj, d, &adof));
4365f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
4375f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]));
4385f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionSetDof(rootSectionAdj, d, adof));
439cb1e1211SMatthew G Knepley     }
440cb1e1211SMatthew G Knepley   }
441cb1e1211SMatthew G Knepley   /* Debugging */
442cb1e1211SMatthew G Knepley   if (debug) {
443cb1e1211SMatthew G Knepley     IS tmp;
4445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n"));
4455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionView(rootSectionAdj, NULL));
4465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Root adjacency indices after compression\n"));
4475f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
4485f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(tmp, NULL));
4495f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&tmp));
450cb1e1211SMatthew G Knepley   }
451cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
4525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd));
4535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionCreate(comm, &sectionAdj));
4545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd));
455cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
45676185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
457cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
458cb1e1211SMatthew G Knepley 
4595f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
4605f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetConstraintDof(section, p, &cdof));
4615f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(section, p, &off));
4625f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(sectionGlobal, p, &goff));
463cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
464cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
465cb1e1211SMatthew G Knepley 
4665f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(leafSectionAdj, off+d, &ldof));
4675f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(rootSectionAdj, off+d, &rdof));
468cb1e1211SMatthew G Knepley       if (ldof > 0) {
469cb1e1211SMatthew G Knepley         /* We do not own this point */
470cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
4715f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionSetDof(sectionAdj, goff+d, rdof));
472cb1e1211SMatthew G Knepley       } else {
473cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
474cb1e1211SMatthew G Knepley       }
475cb1e1211SMatthew G Knepley     }
476cb1e1211SMatthew G Knepley     if (found) continue;
4775f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
4785f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(sectionGlobal, p, &goff));
4795f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
480cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
481cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
482cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof, noff;
483cb1e1211SMatthew G Knepley 
484cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
4855f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(section, padj, &ndof));
4865f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetConstraintDof(section, padj, &ncdof));
4875f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(section, padj, &noff));
488cb1e1211SMatthew G Knepley       for (d = goff; d < goff+dof-cdof; ++d) {
4895f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionAddDof(sectionAdj, d, ndof-ncdof));
490cb1e1211SMatthew G Knepley       }
491cb1e1211SMatthew G Knepley     }
4925f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
49376185916SToby Isaac     if (anDof) {
49476185916SToby Isaac       for (d = goff; d < goff+dof-cdof; ++d) {
4955f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionAddDof(sectionAdj, d, anDof));
49676185916SToby Isaac       }
49776185916SToby Isaac     }
498cb1e1211SMatthew G Knepley   }
4995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionSetUp(sectionAdj));
500cb1e1211SMatthew G Knepley   if (debug) {
5015f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Adjacency Section for Preallocation:\n"));
5025f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionView(sectionAdj, NULL));
503cb1e1211SMatthew G Knepley   }
504cb1e1211SMatthew G Knepley   /* Get adjacent indices */
5055f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetStorageSize(sectionAdj, &numCols));
5065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(numCols, &cols));
507cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
50876185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
509cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
510cb1e1211SMatthew G Knepley 
5115f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(section, p, &dof));
5125f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetConstraintDof(section, p, &cdof));
5135f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(section, p, &off));
5145f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(sectionGlobal, p, &goff));
515cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
516cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
517cb1e1211SMatthew G Knepley 
5185f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(leafSectionAdj, off+d, &ldof));
5195f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(rootSectionAdj, off+d, &rdof));
520cb1e1211SMatthew G Knepley       if (ldof > 0) {
521cb1e1211SMatthew G Knepley         /* We do not own this point */
522cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
523cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
524cb1e1211SMatthew G Knepley 
5255f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(sectionAdj, goff+d, &aoff));
5265f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(rootSectionAdj, off+d, &roff));
5275f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscArraycpy(&cols[aoff], &rootAdj[roff], rdof));
528cb1e1211SMatthew G Knepley       } else {
529cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
530cb1e1211SMatthew G Knepley       }
531cb1e1211SMatthew G Knepley     }
532cb1e1211SMatthew G Knepley     if (found) continue;
5335f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
5345f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
5355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
536cb1e1211SMatthew G Knepley     for (d = goff; d < goff+dof-cdof; ++d) {
537cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
538cb1e1211SMatthew G Knepley 
5395f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(sectionAdj, d, &adof));
5405f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(sectionAdj, d, &aoff));
541cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
542cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
543cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
544cb1e1211SMatthew G Knepley         const PetscInt *ncind;
545cb1e1211SMatthew G Knepley 
546cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
547cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
5485f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetDof(section, padj, &ndof));
5495f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetConstraintDof(section, padj, &ncdof));
5505f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetConstraintIndices(section, padj, &ncind));
5515f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(sectionGlobal, padj, &ngoff));
552cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd, ++i) {
553cb1e1211SMatthew G Knepley           cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
554cb1e1211SMatthew G Knepley         }
555cb1e1211SMatthew G Knepley       }
55676185916SToby Isaac       for (q = 0; q < anDof; q++, i++) {
55776185916SToby Isaac         cols[aoff+i] = anchorAdj[anOff + q];
55876185916SToby Isaac       }
5592c71b3e2SJacob 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);
560cb1e1211SMatthew G Knepley     }
561cb1e1211SMatthew G Knepley   }
5625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&anchorSectionAdj));
5635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&leafSectionAdj));
5645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionDestroy(&rootSectionAdj));
5655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(anchorAdj));
5665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rootAdj));
5675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(tmpAdj));
568cb1e1211SMatthew G Knepley   /* Debugging */
569cb1e1211SMatthew G Knepley   if (debug) {
570cb1e1211SMatthew G Knepley     IS tmp;
5715f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Column indices\n"));
5725f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp));
5735f80ce2aSJacob Faibussowitsch     CHKERRQ(ISView(tmp, NULL));
5745f80ce2aSJacob Faibussowitsch     CHKERRQ(ISDestroy(&tmp));
575cb1e1211SMatthew G Knepley   }
576a9fb6443SMatthew G. Knepley 
577a9fb6443SMatthew G. Knepley   *sA     = sectionAdj;
578a9fb6443SMatthew G. Knepley   *colIdx = cols;
579a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
580a9fb6443SMatthew G. Knepley }
581a9fb6443SMatthew G. Knepley 
58225e402d2SMatthew 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[])
583a9fb6443SMatthew G. Knepley {
584e101f074SMatthew G. Knepley   PetscSection   section;
585a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p;
586a9fb6443SMatthew G. Knepley 
587a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
588a9fb6443SMatthew G. Knepley   /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */
5895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
5902c71b3e2SJacob 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);
591a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
5925f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalSection(dm, &section));
5935f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetChart(section, &pStart, &pEnd));
594a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
595294b7009SMatthew G. Knepley       PetscInt rS, rE;
596a9fb6443SMatthew G. Knepley 
5975f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
598a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
599a9fb6443SMatthew G. Knepley         PetscInt numCols, cStart, c;
600a9fb6443SMatthew G. Knepley 
6015f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetDof(sectionAdj, r, &numCols));
6025f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(sectionAdj, r, &cStart));
603a9fb6443SMatthew G. Knepley         for (c = cStart; c < cStart+numCols; ++c) {
604a9fb6443SMatthew G. Knepley           if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
605a9fb6443SMatthew G. Knepley             ++dnz[r-rStart];
606a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++dnzu[r-rStart];
607a9fb6443SMatthew G. Knepley           } else {
608a9fb6443SMatthew G. Knepley             ++onz[r-rStart];
609a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++onzu[r-rStart];
610a9fb6443SMatthew G. Knepley           }
611a9fb6443SMatthew G. Knepley         }
612a9fb6443SMatthew G. Knepley       }
613a9fb6443SMatthew G. Knepley     }
614a9fb6443SMatthew G. Knepley   } else {
615cb1e1211SMatthew G Knepley     /* Only loop over blocks of rows */
616cb1e1211SMatthew G Knepley     for (r = rStart/bs; r < rEnd/bs; ++r) {
617cb1e1211SMatthew G Knepley       const PetscInt row = r*bs;
618cb1e1211SMatthew G Knepley       PetscInt       numCols, cStart, c;
619cb1e1211SMatthew G Knepley 
6205f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(sectionAdj, row, &numCols));
6215f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(sectionAdj, row, &cStart));
622cb1e1211SMatthew G Knepley       for (c = cStart; c < cStart+numCols; ++c) {
623e7bcfa36SMatthew G. Knepley         if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
624e7bcfa36SMatthew G. Knepley           ++dnz[r-rStart/bs];
625e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++dnzu[r-rStart/bs];
626cb1e1211SMatthew G Knepley         } else {
627e7bcfa36SMatthew G. Knepley           ++onz[r-rStart/bs];
628e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++onzu[r-rStart/bs];
629cb1e1211SMatthew G Knepley         }
630cb1e1211SMatthew G Knepley       }
631cb1e1211SMatthew G Knepley     }
632a9fb6443SMatthew G. Knepley     for (r = 0; r < (rEnd - rStart)/bs; ++r) {
633cb1e1211SMatthew G Knepley       dnz[r]  /= bs;
634cb1e1211SMatthew G Knepley       onz[r]  /= bs;
635cb1e1211SMatthew G Knepley       dnzu[r] /= bs;
636cb1e1211SMatthew G Knepley       onzu[r] /= bs;
637cb1e1211SMatthew G Knepley     }
638cb1e1211SMatthew G Knepley   }
639a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
640a9fb6443SMatthew G. Knepley }
641a9fb6443SMatthew G. Knepley 
64225e402d2SMatthew G. Knepley static PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A)
643a9fb6443SMatthew G. Knepley {
644e101f074SMatthew G. Knepley   PetscSection   section;
645a9fb6443SMatthew G. Knepley   PetscScalar   *values;
646a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0;
647a9fb6443SMatthew G. Knepley 
648a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
6495f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
650a9fb6443SMatthew G. Knepley   for (r = rStart; r < rEnd; ++r) {
6515f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetDof(sectionAdj, r, &len));
652a9fb6443SMatthew G. Knepley     maxRowLen = PetscMax(maxRowLen, len);
653a9fb6443SMatthew G. Knepley   }
6545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscCalloc1(maxRowLen, &values));
655a9fb6443SMatthew G. Knepley   if (f >=0 && bs == 1) {
6565f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalSection(dm, &section));
6575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionGetChart(section, &pStart, &pEnd));
658a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
659294b7009SMatthew G. Knepley       PetscInt rS, rE;
660a9fb6443SMatthew G. Knepley 
6615f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE));
662a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
6637e01ee02SMatthew G. Knepley         PetscInt numCols, cStart;
664a9fb6443SMatthew G. Knepley 
6655f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetDof(sectionAdj, r, &numCols));
6665f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSectionGetOffset(sectionAdj, r, &cStart));
6675f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
668a9fb6443SMatthew G. Knepley       }
669a9fb6443SMatthew G. Knepley     }
670a9fb6443SMatthew G. Knepley   } else {
671a9fb6443SMatthew G. Knepley     for (r = rStart; r < rEnd; ++r) {
672a9fb6443SMatthew G. Knepley       PetscInt numCols, cStart;
673a9fb6443SMatthew G. Knepley 
6745f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetDof(sectionAdj, r, &numCols));
6755f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSectionGetOffset(sectionAdj, r, &cStart));
6765f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
677a9fb6443SMatthew G. Knepley     }
678a9fb6443SMatthew G. Knepley   }
6795f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(values));
680a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
681a9fb6443SMatthew G. Knepley }
682a9fb6443SMatthew G. Knepley 
68325e402d2SMatthew G. Knepley /*@C
68425e402d2SMatthew G. Knepley   DMPlexPreallocateOperator - Calculate the matrix nonzero pattern based upon the information in the DM,
68525e402d2SMatthew G. Knepley   the PetscDS it contains, and the default PetscSection.
68625e402d2SMatthew G. Knepley 
68725e402d2SMatthew G. Knepley   Collective
68825e402d2SMatthew G. Knepley 
6894165533cSJose E. Roman   Input Parameters:
69025e402d2SMatthew G. Knepley + dm   - The DMPlex
69125e402d2SMatthew G. Knepley . bs   - The matrix blocksize
69225e402d2SMatthew G. Knepley . dnz  - An array to hold the number of nonzeros in the diagonal block
69325e402d2SMatthew G. Knepley . onz  - An array to hold the number of nonzeros in the off-diagonal block
69425e402d2SMatthew G. Knepley . dnzu - An array to hold the number of nonzeros in the upper triangle of the diagonal block
69525e402d2SMatthew G. Knepley . onzu - An array to hold the number of nonzeros in the upper triangle of the off-diagonal block
696a2b725a8SWilliam Gropp - fillMatrix - If PETSC_TRUE, fill the matrix with zeros
69725e402d2SMatthew G. Knepley 
6984165533cSJose E. Roman   Output Parameter:
69925e402d2SMatthew G. Knepley . A - The preallocated matrix
70025e402d2SMatthew G. Knepley 
70125e402d2SMatthew G. Knepley   Level: advanced
70225e402d2SMatthew G. Knepley 
70325e402d2SMatthew G. Knepley .seealso: DMCreateMatrix()
70425e402d2SMatthew G. Knepley @*/
705e101f074SMatthew G. Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
706a9fb6443SMatthew G. Knepley {
707a9fb6443SMatthew G. Knepley   MPI_Comm       comm;
708a9fb6443SMatthew G. Knepley   PetscDS        prob;
709a9fb6443SMatthew G. Knepley   MatType        mtype;
710a9fb6443SMatthew G. Knepley   PetscSF        sf, sfDof;
711e101f074SMatthew G. Knepley   PetscSection   section;
712a9fb6443SMatthew G. Knepley   PetscInt      *remoteOffsets;
713a9fb6443SMatthew G. Knepley   PetscSection   sectionAdj[4] = {NULL, NULL, NULL, NULL};
714a9fb6443SMatthew G. Knepley   PetscInt      *cols[4]       = {NULL, NULL, NULL, NULL};
715a9fb6443SMatthew G. Knepley   PetscBool      useCone, useClosure;
7167e01ee02SMatthew G. Knepley   PetscInt       Nf, f, idx, locRows;
717a9fb6443SMatthew G. Knepley   PetscLayout    rLayout;
718a9fb6443SMatthew G. Knepley   PetscBool      isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE;
7196cded2eaSMatthew G. Knepley   PetscMPIInt    size;
720a9fb6443SMatthew G. Knepley 
721a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
722a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
723064a246eSJacob Faibussowitsch   PetscValidHeaderSpecific(A, MAT_CLASSID, 7);
724*dadcf809SJacob Faibussowitsch   if (dnz) PetscValidIntPointer(dnz,3);
725*dadcf809SJacob Faibussowitsch   if (onz) PetscValidIntPointer(onz,4);
726*dadcf809SJacob Faibussowitsch   if (dnzu) PetscValidIntPointer(dnzu,5);
727*dadcf809SJacob Faibussowitsch   if (onzu) PetscValidIntPointer(onzu,6);
7285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDS(dm, &prob));
7295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetPointSF(dm, &sf));
7305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalSection(dm, &section));
7315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL));
7325f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject) dm, &comm));
7335f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(comm, &size));
7345f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0));
735a9fb6443SMatthew G. Knepley   /* Create dof SF based on point SF */
736a9fb6443SMatthew G. Knepley   if (debug) {
737e101f074SMatthew G. Knepley     PetscSection section, sectionGlobal;
738a9fb6443SMatthew G. Knepley     PetscSF      sf;
739a9fb6443SMatthew G. Knepley 
7405f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetPointSF(dm, &sf));
7415f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetLocalSection(dm, &section));
7425f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetGlobalSection(dm, &sectionGlobal));
7435f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Input Section for Preallocation:\n"));
7445f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionView(section, NULL));
7455f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Input Global Section for Preallocation:\n"));
7465f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSectionView(sectionGlobal, NULL));
7476cded2eaSMatthew G. Knepley     if (size > 1) {
7485f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscPrintf(comm, "Input SF for Preallocation:\n"));
7495f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscSFView(sf, NULL));
750a9fb6443SMatthew G. Knepley     }
7516cded2eaSMatthew G. Knepley   }
7525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets));
7535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof));
7545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(remoteOffsets));
7556cded2eaSMatthew G. Knepley   if (debug && size > 1) {
7565f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm, "Dof SF for Preallocation:\n"));
7575f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSFView(sfDof, NULL));
758a9fb6443SMatthew G. Knepley   }
759a9fb6443SMatthew G. Knepley   /* Create allocation vectors from adjacency graph */
7605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetLocalSize(A, &locRows, NULL));
7615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutCreate(comm, &rLayout));
7625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetLocalSize(rLayout, locRows));
7635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetBlockSize(rLayout, 1));
7645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutSetUp(rLayout));
765a9fb6443SMatthew G. Knepley   /* There are 4 types of adjacency */
7665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetNumFields(section, &Nf));
767a9fb6443SMatthew G. Knepley   if (Nf < 1 || bs > 1) {
7685f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGetBasicAdjacency(dm, &useCone, &useClosure));
769a9fb6443SMatthew G. Knepley     idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7705f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7715f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
772a9fb6443SMatthew G. Knepley   } else {
773a9fb6443SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
7745f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetAdjacency(dm, f, &useCone, &useClosure));
775a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7765f80ce2aSJacob Faibussowitsch       if (!sectionAdj[idx]) CHKERRQ(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
7775f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
778a9fb6443SMatthew G. Knepley     }
779a9fb6443SMatthew G. Knepley   }
7805f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFDestroy(&sfDof));
781cb1e1211SMatthew G Knepley   /* Set matrix pattern */
7825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu));
7835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
78489545effSMatthew G. Knepley   /* Check for symmetric storage */
7855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetType(A, &mtype));
7865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
7875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
7885f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
7895f80ce2aSJacob Faibussowitsch   if (isSymBlock || isSymSeqBlock || isSymMPIBlock) CHKERRQ(MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE));
790cb1e1211SMatthew G Knepley   /* Fill matrix with zeros */
791cb1e1211SMatthew G Knepley   if (fillMatrix) {
792a9fb6443SMatthew G. Knepley     if (Nf < 1 || bs > 1) {
7935f80ce2aSJacob Faibussowitsch       CHKERRQ(DMGetBasicAdjacency(dm, &useCone, &useClosure));
794a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
7955f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A));
796a9fb6443SMatthew G. Knepley     } else {
797a9fb6443SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
7985f80ce2aSJacob Faibussowitsch         CHKERRQ(DMGetAdjacency(dm, f, &useCone, &useClosure));
799a9fb6443SMatthew G. Knepley         idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
8005f80ce2aSJacob Faibussowitsch         CHKERRQ(DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A));
801cb1e1211SMatthew G Knepley       }
802cb1e1211SMatthew G Knepley     }
8035f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
8045f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
805cb1e1211SMatthew G Knepley   }
8065f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLayoutDestroy(&rLayout));
8075f80ce2aSJacob Faibussowitsch   for (idx = 0; idx < 4; ++idx) {CHKERRQ(PetscSectionDestroy(&sectionAdj[idx])); CHKERRQ(PetscFree(cols[idx]));}
8085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogEventEnd(DMPLEX_Preallocate,dm,0,0,0));
809cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
810cb1e1211SMatthew G Knepley }
811cb1e1211SMatthew G Knepley 
812cb1e1211SMatthew G Knepley #if 0
813cb1e1211SMatthew 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)
814cb1e1211SMatthew G Knepley {
815cb1e1211SMatthew G Knepley   PetscInt       *tmpClosure,*tmpAdj,*visits;
816cb1e1211SMatthew G Knepley   PetscInt        c,cStart,cEnd,pStart,pEnd;
817cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
818cb1e1211SMatthew G Knepley 
819cb1e1211SMatthew G Knepley   PetscFunctionBegin;
8205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetDimension(dm, &dim));
8215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepth(dm, &depth));
8225f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
8265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSectionGetChart(section, &pStart, &pEnd));
827cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
828cb1e1211SMatthew G Knepley 
8295f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits));
8305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(lvisits,pEnd-pStart));
8315f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArrayzero(visits,pEnd-pStart));
8325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
833cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
834cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
8355f80ce2aSJacob Faibussowitsch     CHKERRQ(DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support));
836cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
837cb1e1211SMatthew G Knepley   }
8385f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM));
8395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM));
8405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE));
8415f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits));
842cb1e1211SMatthew G Knepley 
8435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFGetRootRanks());
844cb1e1211SMatthew G Knepley 
8455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner));
846cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
8475f80ce2aSJacob Faibussowitsch     CHKERRQ(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 
8555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM));
8565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM));
857cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
858cb1e1211SMatthew G Knepley }
859cb1e1211SMatthew G Knepley #endif
860