xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision 9566063d113dddea24716c546802770db7481bc0)
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;
14*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
15*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
16*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) section), &adjSec));
17*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section,&pStart,&pEnd));
18*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(adjSec,pStart,pEnd));
1976185916SToby Isaac 
20*9566063dSJacob 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 */
28*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)aSec),&inverseSec));
29*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(inverseSec,pStart,pEnd));
30*9566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(aIS, &aSize));
31*9566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(aIS, &anchors));
3276185916SToby Isaac 
3376185916SToby Isaac     for (p = 0; p < aSize; p++) {
3476185916SToby Isaac       PetscInt a = anchors[p];
3576185916SToby Isaac 
36*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionAddDof(inverseSec,a,1));
3776185916SToby Isaac     }
38*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(inverseSec));
39*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(inverseSec,&iSize));
40*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(iSize,&inverse));
41*9566063dSJacob Faibussowitsch     PetscCall(PetscCalloc1(pEnd-pStart,&offsets));
42*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(aSec,&aStart,&aEnd));
4376185916SToby Isaac     for (p = aStart; p < aEnd; p++) {
4476185916SToby Isaac       PetscInt dof, off;
4576185916SToby Isaac 
46*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(aSec, p, &dof));
47*9566063dSJacob 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];
53*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(inverseSec, a, &iOff));
5476185916SToby Isaac         inverse[iOff + offsets[a-pStart]++] = p;
5576185916SToby Isaac       }
5676185916SToby Isaac     }
57*9566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(aIS, &anchors));
58*9566063dSJacob 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 
80*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(inverseSec,p,&iDof));
8176185916SToby Isaac       if (!iDof) continue;
82*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(inverseSec,p,&iOff));
83*9566063dSJacob 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];
88*9566063dSJacob 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;
96*9566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section,qAdj,&qAdjDof));
97*9566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof));
9876185916SToby Isaac           iNew += qAdjDof - qAdjCDof;
9976185916SToby Isaac         }
100*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(adjSec,p,iNew));
10176185916SToby Isaac       }
10276185916SToby Isaac     }
10376185916SToby Isaac 
104*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetUp(adjSec));
105*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(adjSec,&adjSize));
106*9566063dSJacob 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 
111*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(inverseSec,p,&iDof));
11276185916SToby Isaac       if (!iDof) continue;
113*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(inverseSec,p,&iOff));
114*9566063dSJacob Faibussowitsch       PetscCall(DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP));
115*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(adjSec,p,&aDof));
116*9566063dSJacob 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];
122*9566063dSJacob 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;
130*9566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetDof(section,qAdj,&qAdjDof));
131*9566063dSJacob Faibussowitsch           PetscCall(PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof));
132*9566063dSJacob 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       }
138*9566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]));
139*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(adjSec,p,aDof));
14076185916SToby Isaac     }
14176185916SToby Isaac     *anchorAdj = adj;
14276185916SToby Isaac 
14376185916SToby Isaac     /* clean up */
144*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&inverseSec));
145*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(inverse));
146*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjP));
147*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(tmpAdjQ));
14876185916SToby Isaac   }
14976185916SToby Isaac   else {
15076185916SToby Isaac     *anchorAdj = NULL;
151*9566063dSJacob 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;
172*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
173*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL));
174*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
175*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
176*9566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
177*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
178*9566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
179*9566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL));
18047a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
181*9566063dSJacob Faibussowitsch   PetscCallMPI(MPIU_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm));
182cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
183*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
184*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &numDof));
185*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &leafSectionAdj));
186*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(leafSectionAdj, 0, numDof));
187*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &rootSectionAdj));
188*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(rootSectionAdj, 0, numDof));
189cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
190*9566063dSJacob 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   */
218*9566063dSJacob 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;
224*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
225*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
226*9566063dSJacob 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;
232*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
233*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
234cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
235*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof));
236cb1e1211SMatthew G Knepley       }
237cb1e1211SMatthew G Knepley     }
238*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
23976185916SToby Isaac     if (anDof) {
24076185916SToby Isaac       for (d = off; d < off+dof; ++d) {
241*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(leafSectionAdj, d, anDof));
24276185916SToby Isaac       }
24376185916SToby Isaac     }
244cb1e1211SMatthew G Knepley   }
245*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
246cb1e1211SMatthew G Knepley   if (debug) {
247*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n"));
248*9566063dSJacob 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) {
252*9566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
253*9566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM));
254cb1e1211SMatthew G Knepley   }
255cb1e1211SMatthew G Knepley   if (debug) {
256*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n"));
257*9566063dSJacob Faibussowitsch     PetscCall(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 
263*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
264*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
265cb1e1211SMatthew G Knepley     if (!dof) continue;
266*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
267cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
268*9566063dSJacob Faibussowitsch     PetscCall(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;
274*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
275*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
276cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
277*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof));
278cb1e1211SMatthew G Knepley       }
279cb1e1211SMatthew G Knepley     }
280*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
28176185916SToby Isaac     if (anDof) {
28276185916SToby Isaac       for (d = off; d < off+dof; ++d) {
283*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(rootSectionAdj, d, anDof));
28476185916SToby Isaac       }
28576185916SToby Isaac     }
286cb1e1211SMatthew G Knepley   }
287*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSectionAdj));
288cb1e1211SMatthew G Knepley   if (debug) {
289*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n"));
290*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
291cb1e1211SMatthew G Knepley   }
292cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
293*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets));
294*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj));
295*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
2966cded2eaSMatthew G. Knepley   if (debug && size > 1) {
297*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency SF for Preallocation:\n"));
298*9566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfAdj, NULL));
299cb1e1211SMatthew G Knepley   }
300cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
301*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(leafSectionAdj));
302*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(leafSectionAdj, &adjSize));
303*9566063dSJacob Faibussowitsch   PetscCall(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;
309*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
310*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
311*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
312*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
313*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
314cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
315cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
316cb1e1211SMatthew G Knepley 
317*9566063dSJacob Faibussowitsch       PetscCall(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;
323*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
324*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
325*9566063dSJacob Faibussowitsch         PetscCall(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;
340*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Leaf adjacency indices\n"));
341*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp));
342*9566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
343*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
344cb1e1211SMatthew G Knepley   }
345543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
346*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSectionAdj, &adjSize));
347*9566063dSJacob Faibussowitsch   PetscCall(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 
353*9566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeBegin(sfAdj, &indegree));
354*9566063dSJacob Faibussowitsch     PetscCall(PetscSFComputeDegreeEnd(sfAdj, &indegree));
355543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
356*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(radjsize, &remoteadj));
357*9566063dSJacob Faibussowitsch     PetscCall(PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj));
358*9566063dSJacob Faibussowitsch     PetscCall(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);
365*9566063dSJacob Faibussowitsch     PetscCall(PetscFree(remoteadj));
366cb1e1211SMatthew G Knepley   }
367*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfAdj));
368*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
369cb1e1211SMatthew G Knepley   /* Debugging */
370cb1e1211SMatthew G Knepley   if (debug) {
371cb1e1211SMatthew G Knepley     IS tmp;
372*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after gather\n"));
373*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
374*9566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
375*9566063dSJacob Faibussowitsch     PetscCall(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 
381*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
382*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
383cb1e1211SMatthew G Knepley     if (!dof) continue;
384*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
385cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
386*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
387*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
388*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff));
389cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
390cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
391cb1e1211SMatthew G Knepley 
392*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
393*9566063dSJacob Faibussowitsch       PetscCall(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;
404*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
405*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
406*9566063dSJacob Faibussowitsch         PetscCall(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;
417*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices\n"));
418*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
419*9566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
420*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
421cb1e1211SMatthew G Knepley   }
422cb1e1211SMatthew G Knepley   /* Compress indices */
423*9566063dSJacob Faibussowitsch   PetscCall(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 
428*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
429*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
430*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
431cb1e1211SMatthew G Knepley     if (!dof) continue;
432*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof));
433cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
434cb1e1211SMatthew G Knepley     for (d = off; d < off+dof-cdof; ++d) {
435*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof));
436*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff));
437*9566063dSJacob Faibussowitsch       PetscCall(PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]));
438*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionSetDof(rootSectionAdj, d, adof));
439cb1e1211SMatthew G Knepley     }
440cb1e1211SMatthew G Knepley   }
441cb1e1211SMatthew G Knepley   /* Debugging */
442cb1e1211SMatthew G Knepley   if (debug) {
443cb1e1211SMatthew G Knepley     IS tmp;
444*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n"));
445*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(rootSectionAdj, NULL));
446*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Root adjacency indices after compression\n"));
447*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp));
448*9566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
449*9566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&tmp));
450cb1e1211SMatthew G Knepley   }
451cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
452*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd));
453*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &sectionAdj));
454*9566063dSJacob Faibussowitsch   PetscCall(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 
459*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
460*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
461*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
462*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
463cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
464cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
465cb1e1211SMatthew G Knepley 
466*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off+d, &ldof));
467*9566063dSJacob Faibussowitsch       PetscCall(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) {
471*9566063dSJacob Faibussowitsch         PetscCall(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;
477*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
478*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
479*9566063dSJacob Faibussowitsch     PetscCall(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;
485*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(section, padj, &ndof));
486*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
487*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(section, padj, &noff));
488cb1e1211SMatthew G Knepley       for (d = goff; d < goff+dof-cdof; ++d) {
489*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(sectionAdj, d, ndof-ncdof));
490cb1e1211SMatthew G Knepley       }
491cb1e1211SMatthew G Knepley     }
492*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
49376185916SToby Isaac     if (anDof) {
49476185916SToby Isaac       for (d = goff; d < goff+dof-cdof; ++d) {
495*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionAddDof(sectionAdj, d, anDof));
49676185916SToby Isaac       }
49776185916SToby Isaac     }
498cb1e1211SMatthew G Knepley   }
499*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(sectionAdj));
500cb1e1211SMatthew G Knepley   if (debug) {
501*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation:\n"));
502*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionAdj, NULL));
503cb1e1211SMatthew G Knepley   }
504cb1e1211SMatthew G Knepley   /* Get adjacent indices */
505*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(sectionAdj, &numCols));
506*9566063dSJacob Faibussowitsch   PetscCall(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 
511*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
512*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(section, p, &cdof));
513*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(section, p, &off));
514*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff));
515cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
516cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
517cb1e1211SMatthew G Knepley 
518*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSectionAdj, off+d, &ldof));
519*9566063dSJacob Faibussowitsch       PetscCall(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 
525*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, goff+d, &aoff));
526*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(rootSectionAdj, off+d, &roff));
527*9566063dSJacob Faibussowitsch         PetscCall(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;
533*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj));
534*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof));
535*9566063dSJacob Faibussowitsch     PetscCall(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 
539*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, d, &adof));
540*9566063dSJacob Faibussowitsch       PetscCall(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;
548*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(section, padj, &ndof));
549*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof));
550*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintIndices(section, padj, &ncind));
551*9566063dSJacob Faibussowitsch         PetscCall(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   }
562*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&anchorSectionAdj));
563*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&leafSectionAdj));
564*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&rootSectionAdj));
565*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(anchorAdj));
566*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(rootAdj));
567*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(tmpAdj));
568cb1e1211SMatthew G Knepley   /* Debugging */
569cb1e1211SMatthew G Knepley   if (debug) {
570cb1e1211SMatthew G Knepley     IS tmp;
571*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Column indices\n"));
572*9566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp));
573*9566063dSJacob Faibussowitsch     PetscCall(ISView(tmp, NULL));
574*9566063dSJacob Faibussowitsch     PetscCall(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 */
589*9566063dSJacob Faibussowitsch   PetscCall(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) {
592*9566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
593*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
594a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
595294b7009SMatthew G. Knepley       PetscInt rS, rE;
596a9fb6443SMatthew G. Knepley 
597*9566063dSJacob Faibussowitsch       PetscCall(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 
601*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
602*9566063dSJacob Faibussowitsch         PetscCall(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 
620*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, row, &numCols));
621*9566063dSJacob Faibussowitsch       PetscCall(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;
649*9566063dSJacob Faibussowitsch   PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd));
650a9fb6443SMatthew G. Knepley   for (r = rStart; r < rEnd; ++r) {
651*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(sectionAdj, r, &len));
652a9fb6443SMatthew G. Knepley     maxRowLen = PetscMax(maxRowLen, len);
653a9fb6443SMatthew G. Knepley   }
654*9566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(maxRowLen, &values));
655a9fb6443SMatthew G. Knepley   if (f >=0 && bs == 1) {
656*9566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
657*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
658a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
659294b7009SMatthew G. Knepley       PetscInt rS, rE;
660a9fb6443SMatthew G. Knepley 
661*9566063dSJacob Faibussowitsch       PetscCall(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 
665*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
666*9566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
667*9566063dSJacob Faibussowitsch         PetscCall(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 
674*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols));
675*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart));
676*9566063dSJacob Faibussowitsch       PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES));
677a9fb6443SMatthew G. Knepley     }
678a9fb6443SMatthew G. Knepley   }
679*9566063dSJacob Faibussowitsch   PetscCall(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);
724dadcf809SJacob Faibussowitsch   if (dnz) PetscValidIntPointer(dnz,3);
725dadcf809SJacob Faibussowitsch   if (onz) PetscValidIntPointer(onz,4);
726dadcf809SJacob Faibussowitsch   if (dnzu) PetscValidIntPointer(dnzu,5);
727dadcf809SJacob Faibussowitsch   if (onzu) PetscValidIntPointer(onzu,6);
728*9566063dSJacob Faibussowitsch   PetscCall(DMGetDS(dm, &prob));
729*9566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sf));
730*9566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
731*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL));
732*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
733*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
734*9566063dSJacob Faibussowitsch   PetscCall(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 
740*9566063dSJacob Faibussowitsch     PetscCall(DMGetPointSF(dm, &sf));
741*9566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
742*9566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
743*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Section for Preallocation:\n"));
744*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(section, NULL));
745*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Input Global Section for Preallocation:\n"));
746*9566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(sectionGlobal, NULL));
7476cded2eaSMatthew G. Knepley     if (size > 1) {
748*9566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Input SF for Preallocation:\n"));
749*9566063dSJacob Faibussowitsch       PetscCall(PetscSFView(sf, NULL));
750a9fb6443SMatthew G. Knepley     }
7516cded2eaSMatthew G. Knepley   }
752*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets));
753*9566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof));
754*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
7556cded2eaSMatthew G. Knepley   if (debug && size > 1) {
756*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Dof SF for Preallocation:\n"));
757*9566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfDof, NULL));
758a9fb6443SMatthew G. Knepley   }
759a9fb6443SMatthew G. Knepley   /* Create allocation vectors from adjacency graph */
760*9566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A, &locRows, NULL));
761*9566063dSJacob Faibussowitsch   PetscCall(PetscLayoutCreate(comm, &rLayout));
762*9566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetLocalSize(rLayout, locRows));
763*9566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetBlockSize(rLayout, 1));
764*9566063dSJacob Faibussowitsch   PetscCall(PetscLayoutSetUp(rLayout));
765a9fb6443SMatthew G. Knepley   /* There are 4 types of adjacency */
766*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetNumFields(section, &Nf));
767a9fb6443SMatthew G. Knepley   if (Nf < 1 || bs > 1) {
768*9566063dSJacob Faibussowitsch     PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
769a9fb6443SMatthew G. Knepley     idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
770*9566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
771*9566063dSJacob Faibussowitsch     PetscCall(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) {
774*9566063dSJacob Faibussowitsch       PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
775a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
776*9566063dSJacob Faibussowitsch       if (!sectionAdj[idx]) PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]));
777*9566063dSJacob Faibussowitsch       PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu));
778a9fb6443SMatthew G. Knepley     }
779a9fb6443SMatthew G. Knepley   }
780*9566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfDof));
781cb1e1211SMatthew G Knepley   /* Set matrix pattern */
782*9566063dSJacob Faibussowitsch   PetscCall(MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu));
783*9566063dSJacob Faibussowitsch   PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
78489545effSMatthew G. Knepley   /* Check for symmetric storage */
785*9566063dSJacob Faibussowitsch   PetscCall(MatGetType(A, &mtype));
786*9566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock));
787*9566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock));
788*9566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock));
789*9566063dSJacob Faibussowitsch   if (isSymBlock || isSymSeqBlock || isSymMPIBlock) PetscCall(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) {
793*9566063dSJacob Faibussowitsch       PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
794a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
795*9566063dSJacob Faibussowitsch       PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A));
796a9fb6443SMatthew G. Knepley     } else {
797a9fb6443SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
798*9566063dSJacob Faibussowitsch         PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure));
799a9fb6443SMatthew G. Knepley         idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
800*9566063dSJacob Faibussowitsch         PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A));
801cb1e1211SMatthew G Knepley       }
802cb1e1211SMatthew G Knepley     }
803*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
804*9566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
805cb1e1211SMatthew G Knepley   }
806*9566063dSJacob Faibussowitsch   PetscCall(PetscLayoutDestroy(&rLayout));
807*9566063dSJacob Faibussowitsch   for (idx = 0; idx < 4; ++idx) {PetscCall(PetscSectionDestroy(&sectionAdj[idx])); PetscCall(PetscFree(cols[idx]));}
808*9566063dSJacob Faibussowitsch   PetscCall(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 
818cb1e1211SMatthew G Knepley   PetscFunctionBegin;
819*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
820*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
821*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize));
822cb1e1211SMatthew G Knepley 
823e7b80042SMatthew G. Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1));
824cb1e1211SMatthew G Knepley 
825*9566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
826cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
827cb1e1211SMatthew G Knepley 
828*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits));
829*9566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(lvisits,pEnd-pStart));
830*9566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(visits,pEnd-pStart));
831*9566063dSJacob Faibussowitsch   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
832cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
833cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
834*9566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support));
835cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
836cb1e1211SMatthew G Knepley   }
837*9566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM));
838*9566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM));
839*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE));
840*9566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits));
841cb1e1211SMatthew G Knepley 
842*9566063dSJacob Faibussowitsch   PetscCall(PetscSFGetRootRanks());
843cb1e1211SMatthew G Knepley 
844*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner));
845cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
846*9566063dSJacob Faibussowitsch     PetscCall(PetscArrayzero(cellmat,maxClosureSize*maxClosureSize));
847cb1e1211SMatthew G Knepley     /*
848cb1e1211SMatthew G Knepley      Depth-first walk of transitive closure.
849cb1e1211SMatthew 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.
850cb1e1211SMatthew G Knepley      This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
851cb1e1211SMatthew G Knepley      */
852cb1e1211SMatthew G Knepley   }
853cb1e1211SMatthew G Knepley 
854*9566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM));
855*9566063dSJacob Faibussowitsch   PetscCall(PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM));
856cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
857cb1e1211SMatthew G Knepley }
858cb1e1211SMatthew G Knepley #endif
859