xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision 543482b8f387d3b221f55b823e6773ff9fd99b03)
1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2cb1e1211SMatthew G Knepley #include <petsc-private/isimpl.h>
30c312b8eSJed Brown #include <petscsf.h>
4cb1e1211SMatthew G Knepley 
5cb1e1211SMatthew G Knepley #undef __FUNCT__
676185916SToby Isaac #define __FUNCT__ "DMPlexComputeAnchorAdjacencies"
776185916SToby Isaac /* get adjacencies due to point-to-point constraints that can't be found with DMPlexGetAdjacency() */
876185916SToby Isaac static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscSection section, PetscSection sectionGlobal, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[])
976185916SToby Isaac {
1076185916SToby Isaac   PetscInt       pStart, pEnd;
1176185916SToby Isaac   PetscSection   adjSec, aSec;
1276185916SToby Isaac   IS             aIS;
1376185916SToby Isaac   PetscErrorCode ierr;
1476185916SToby Isaac 
1576185916SToby Isaac   PetscFunctionBegin;
1676185916SToby Isaac 
1776185916SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&adjSec);CHKERRQ(ierr);
1876185916SToby Isaac   ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr);
1976185916SToby Isaac   ierr = PetscSectionSetChart(adjSec,pStart,pEnd);CHKERRQ(ierr);
2076185916SToby Isaac 
21a17985deSToby Isaac   ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
2276185916SToby Isaac   if (aSec) {
2376185916SToby Isaac     const PetscInt *anchors;
2476185916SToby Isaac     PetscInt       p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize;
2576185916SToby Isaac     PetscInt       *tmpAdjP = NULL, *tmpAdjQ = NULL;
2676185916SToby Isaac     PetscSection   inverseSec;
2776185916SToby Isaac 
2876185916SToby Isaac     /* invert the constraint-to-anchor map */
2976185916SToby Isaac     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)aSec),&inverseSec);CHKERRQ(ierr);
3076185916SToby Isaac     ierr = PetscSectionSetChart(inverseSec,pStart,pEnd);CHKERRQ(ierr);
3176185916SToby Isaac     ierr = ISGetLocalSize(aIS, &aSize);CHKERRQ(ierr);
3276185916SToby Isaac     ierr = ISGetIndices(aIS, &anchors);CHKERRQ(ierr);
3376185916SToby Isaac 
3476185916SToby Isaac     for (p = 0; p < aSize; p++) {
3576185916SToby Isaac       PetscInt a = anchors[p];
3676185916SToby Isaac 
3776185916SToby Isaac       ierr = PetscSectionAddDof(inverseSec,a,1);CHKERRQ(ierr);
3876185916SToby Isaac     }
3976185916SToby Isaac     ierr = PetscSectionSetUp(inverseSec);CHKERRQ(ierr);
4076185916SToby Isaac     ierr = PetscSectionGetStorageSize(inverseSec,&iSize);CHKERRQ(ierr);
4176185916SToby Isaac     ierr = PetscMalloc1(iSize,&inverse);CHKERRQ(ierr);
4276185916SToby Isaac     ierr = PetscCalloc1(pEnd-pStart,&offsets);CHKERRQ(ierr);
4376185916SToby Isaac     ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
4476185916SToby Isaac     for (p = aStart; p < aEnd; p++) {
4576185916SToby Isaac       PetscInt dof, off;
4676185916SToby Isaac 
4776185916SToby Isaac       ierr = PetscSectionGetDof(aSec, p, &dof);CHKERRQ(ierr);
4876185916SToby Isaac       ierr = PetscSectionGetOffset(aSec, p, &off);CHKERRQ(ierr);
4976185916SToby Isaac 
5076185916SToby Isaac       for (q = 0; q < dof; q++) {
5176185916SToby Isaac         PetscInt iOff;
5276185916SToby Isaac 
5376185916SToby Isaac         a = anchors[off + q];
5476185916SToby Isaac         ierr = PetscSectionGetOffset(inverseSec, a, &iOff);CHKERRQ(ierr);
5576185916SToby Isaac         inverse[iOff + offsets[a-pStart]++] = p;
5676185916SToby Isaac       }
5776185916SToby Isaac     }
5876185916SToby Isaac     ierr = ISRestoreIndices(aIS, &anchors);CHKERRQ(ierr);
5976185916SToby Isaac     ierr = PetscFree(offsets);CHKERRQ(ierr);
6076185916SToby Isaac 
6176185916SToby Isaac     /* construct anchorAdj and adjSec
6276185916SToby Isaac      *
6376185916SToby Isaac      * loop over anchors:
6476185916SToby Isaac      *   construct anchor adjacency
6576185916SToby Isaac      *   loop over constrained:
6676185916SToby Isaac      *     construct constrained adjacency
6776185916SToby Isaac      *     if not in anchor adjacency, add to dofs
6876185916SToby Isaac      * setup adjSec, allocate anchorAdj
6976185916SToby Isaac      * loop over anchors:
7076185916SToby Isaac      *   construct anchor adjacency
7176185916SToby Isaac      *   loop over constrained:
7276185916SToby Isaac      *     construct constrained adjacency
7376185916SToby Isaac      *     if not in anchor adjacency
7476185916SToby Isaac      *       if not already in list, put in list
7576185916SToby Isaac      *   sort, unique, reduce dof count
7676185916SToby Isaac      * optional: compactify
7776185916SToby Isaac      */
7876185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
7976185916SToby Isaac       PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE;
8076185916SToby Isaac 
8176185916SToby Isaac       ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr);
8276185916SToby Isaac       if (!iDof) continue;
8376185916SToby Isaac       ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr);
8476185916SToby Isaac       ierr = DMPlexGetAdjacency(dm,p,&numAdjP,&tmpAdjP);CHKERRQ(ierr);
8576185916SToby Isaac       for (i = 0; i < iDof; i++) {
8676185916SToby Isaac         PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE;
8776185916SToby Isaac 
8876185916SToby Isaac         q = inverse[iOff + i];
8976185916SToby Isaac         ierr = DMPlexGetAdjacency(dm,q,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr);
9076185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
9176185916SToby Isaac           qAdj = tmpAdjQ[r];
9276185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
9376185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
9476185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
9576185916SToby Isaac           }
9676185916SToby Isaac           if (s < numAdjP) continue;
9776185916SToby Isaac           ierr  = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr);
9876185916SToby Isaac           ierr  = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr);
9976185916SToby Isaac           iNew += qAdjDof - qAdjCDof;
10076185916SToby Isaac         }
10176185916SToby Isaac         ierr = PetscSectionAddDof(adjSec,p,iNew);CHKERRQ(ierr);
10276185916SToby Isaac       }
10376185916SToby Isaac     }
10476185916SToby Isaac 
10576185916SToby Isaac     ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr);
10676185916SToby Isaac     ierr = PetscSectionGetStorageSize(adjSec,&adjSize);CHKERRQ(ierr);
10776185916SToby Isaac     ierr = PetscMalloc1(adjSize,&adj);CHKERRQ(ierr);
10876185916SToby Isaac 
10976185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
11076185916SToby Isaac       PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE;
11176185916SToby Isaac 
11276185916SToby Isaac       ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr);
11376185916SToby Isaac       if (!iDof) continue;
11476185916SToby Isaac       ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr);
11576185916SToby Isaac       ierr = DMPlexGetAdjacency(dm,p,&numAdjP,&tmpAdjP);CHKERRQ(ierr);
11676185916SToby Isaac       ierr = PetscSectionGetDof(adjSec,p,&aDof);CHKERRQ(ierr);
11776185916SToby Isaac       ierr = PetscSectionGetOffset(adjSec,p,&aOff);CHKERRQ(ierr);
11876185916SToby Isaac       aOffOrig = aOff;
11976185916SToby Isaac       for (i = 0; i < iDof; i++) {
12076185916SToby Isaac         PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE;
12176185916SToby Isaac 
12276185916SToby Isaac         q = inverse[iOff + i];
12376185916SToby Isaac         ierr = DMPlexGetAdjacency(dm,q,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr);
12476185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
12576185916SToby Isaac           qAdj = tmpAdjQ[r];
12676185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
12776185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
12876185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
12976185916SToby Isaac           }
13076185916SToby Isaac           if (s < numAdjP) continue;
13176185916SToby Isaac           ierr  = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr);
13276185916SToby Isaac           ierr  = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr);
13376185916SToby Isaac           ierr  = PetscSectionGetOffset(sectionGlobal,qAdj,&qAdjOff);CHKERRQ(ierr);
13476185916SToby Isaac           for (nd = 0; nd < qAdjDof-qAdjCDof; ++nd) {
13576185916SToby Isaac             adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff+1) : qAdjOff) + nd;
13676185916SToby Isaac           }
13776185916SToby Isaac         }
13876185916SToby Isaac       }
13976185916SToby Isaac       ierr = PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]);CHKERRQ(ierr);
14076185916SToby Isaac       ierr = PetscSectionSetDof(adjSec,p,aDof);CHKERRQ(ierr);
14176185916SToby Isaac     }
14276185916SToby Isaac     *anchorAdj = adj;
14376185916SToby Isaac 
14476185916SToby Isaac     /* clean up */
14576185916SToby Isaac     ierr = PetscSectionDestroy(&inverseSec);CHKERRQ(ierr);
14676185916SToby Isaac     ierr = PetscFree(inverse);CHKERRQ(ierr);
14776185916SToby Isaac     ierr = PetscFree(tmpAdjP);CHKERRQ(ierr);
14876185916SToby Isaac     ierr = PetscFree(tmpAdjQ);CHKERRQ(ierr);
14976185916SToby Isaac   }
15076185916SToby Isaac   else {
15176185916SToby Isaac     *anchorAdj = NULL;
15276185916SToby Isaac     ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr);
15376185916SToby Isaac   }
15476185916SToby Isaac   *anchorSectionAdj = adjSec;
15576185916SToby Isaac   PetscFunctionReturn(0);
15676185916SToby Isaac }
15776185916SToby Isaac 
15876185916SToby Isaac #undef __FUNCT__
159cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator"
160cb1e1211SMatthew G Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
161cb1e1211SMatthew G Knepley {
162cb1e1211SMatthew G Knepley   MPI_Comm           comm;
16389545effSMatthew G. Knepley   MatType            mtype;
164cb1e1211SMatthew G Knepley   PetscSF            sf, sfDof, sfAdj;
16576185916SToby Isaac   PetscSection       leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj;
16647a6104aSMatthew G. Knepley   PetscInt           nroots, nleaves, l, p;
167cb1e1211SMatthew G Knepley   const PetscInt    *leaves;
168cb1e1211SMatthew G Knepley   const PetscSFNode *remotes;
169cb1e1211SMatthew G Knepley   PetscInt           dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols;
1708821704fSMatthew G. Knepley   PetscInt          *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets;
17170034214SMatthew G. Knepley   PetscInt           adjSize;
172cb1e1211SMatthew G Knepley   PetscLayout        rLayout;
173cb1e1211SMatthew G Knepley   PetscInt           locRows, rStart, rEnd, r;
174cb1e1211SMatthew G Knepley   PetscMPIInt        size;
17570034214SMatthew G. Knepley   PetscBool          doCommLocal, doComm, debug = PETSC_FALSE, isSymBlock, isSymSeqBlock, isSymMPIBlock;
176e228b242SToby Isaac   PetscBool          useAnchors;
177cb1e1211SMatthew G Knepley   PetscErrorCode     ierr;
178cb1e1211SMatthew G Knepley 
179cb1e1211SMatthew G Knepley   PetscFunctionBegin;
18070034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
18170034214SMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
18270034214SMatthew G. Knepley   PetscValidHeaderSpecific(sectionGlobal, PETSC_SECTION_CLASSID, 4);
18370034214SMatthew G. Knepley   PetscValidHeaderSpecific(A, MAT_CLASSID, 9);
1848c598124SMatthew G. Knepley   if (dnz)  PetscValidPointer(dnz,5);
1858c598124SMatthew G. Knepley   if (onz)  PetscValidPointer(onz,6);
1868c598124SMatthew G. Knepley   if (dnzu) PetscValidPointer(dnzu,7);
1878c598124SMatthew G. Knepley   if (onzu) PetscValidPointer(onzu,8);
188a72f3261SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr);
189cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
190cb1e1211SMatthew G Knepley   ierr = PetscOptionsGetBool(NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr);
191cb1e1211SMatthew G Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
192c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
193cb1e1211SMatthew G Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
19447a6104aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
19547a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
19647a6104aSMatthew G. Knepley   ierr = MPI_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr);
197cb1e1211SMatthew G Knepley   /* Create dof SF based on point SF */
198cb1e1211SMatthew G Knepley   if (debug) {
199cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Input Section for Preallocation:\n");CHKERRQ(ierr);
200cb1e1211SMatthew G Knepley     ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
201cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Input Global Section for Preallocation:\n");CHKERRQ(ierr);
202cb1e1211SMatthew G Knepley     ierr = PetscSectionView(sectionGlobal, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
203cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Input SF for Preallocation:\n");CHKERRQ(ierr);
204cb1e1211SMatthew G Knepley     ierr = PetscSFView(sf, NULL);CHKERRQ(ierr);
205cb1e1211SMatthew G Knepley   }
206cb1e1211SMatthew G Knepley   ierr = PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets);CHKERRQ(ierr);
207cb1e1211SMatthew G Knepley   ierr = PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof);CHKERRQ(ierr);
208cb1e1211SMatthew G Knepley   if (debug) {
209cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Dof SF for Preallocation:\n");CHKERRQ(ierr);
210cb1e1211SMatthew G Knepley     ierr = PetscSFView(sfDof, NULL);CHKERRQ(ierr);
211cb1e1211SMatthew G Knepley   }
212cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
213cb1e1211SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
214cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &numDof);CHKERRQ(ierr);
215cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &leafSectionAdj);CHKERRQ(ierr);
216cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(leafSectionAdj, 0, numDof);CHKERRQ(ierr);
217cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &rootSectionAdj);CHKERRQ(ierr);
218cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(rootSectionAdj, 0, numDof);CHKERRQ(ierr);
219cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
220cb1e1211SMatthew G Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes);CHKERRQ(ierr);
2218b0b4c70SToby Isaac   /* use constraints in finding adjacency in this routine */
222e228b242SToby Isaac   ierr = DMPlexGetAdjacencyUseAnchors(dm,&useAnchors);CHKERRQ(ierr);
223a17985deSToby Isaac   ierr = DMPlexSetAdjacencyUseAnchors(dm,PETSC_TRUE);CHKERRQ(ierr);
224cb1e1211SMatthew G Knepley 
225cb1e1211SMatthew G Knepley   /*
226964bf7afSMatthew G. Knepley    section        - maps points to (# dofs, local dofs)
227964bf7afSMatthew G. Knepley    sectionGlobal  - maps points to (# dofs, global dofs)
228964bf7afSMatthew G. Knepley    leafSectionAdj - maps unowned local dofs to # adj dofs
229964bf7afSMatthew G. Knepley    rootSectionAdj - maps   owned local dofs to # adj dofs
230964bf7afSMatthew G. Knepley    adj            - adj global dofs indexed by leafSectionAdj
231964bf7afSMatthew G. Knepley    rootAdj        - adj global dofs indexed by rootSectionAdj
232964bf7afSMatthew G. Knepley    sf    - describes shared points across procs
233964bf7afSMatthew G. Knepley    sfDof - describes shared dofs across procs
234964bf7afSMatthew G. Knepley    sfAdj - describes shared adjacent dofs across procs
235cb1e1211SMatthew G Knepley    ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
23676185916SToby Isaac   (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj
23776185916SToby Isaac        (This is done in DMPlexComputeAnchorAdjacencies())
238cb1e1211SMatthew G Knepley     1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
239cb1e1211SMatthew G Knepley        Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
240cb1e1211SMatthew G Knepley     2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
241cb1e1211SMatthew G Knepley        Create sfAdj connecting rootSectionAdj and leafSectionAdj
242cb1e1211SMatthew G Knepley     3. Visit unowned points on interface, write adjacencies to adj
243cb1e1211SMatthew G Knepley        Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
244cb1e1211SMatthew G Knepley     4. Visit owned points on interface, write adjacencies to rootAdj
245cb1e1211SMatthew G Knepley        Remove redundancy in rootAdj
246cb1e1211SMatthew G Knepley    ** The last two traversals use transitive closure
247cb1e1211SMatthew G Knepley     5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
248cb1e1211SMatthew G Knepley        Allocate memory addressed by sectionAdj (cols)
249cb1e1211SMatthew G Knepley     6. Visit all owned points in the subdomain, insert dof adjacencies into cols
250cb1e1211SMatthew G Knepley    ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
251cb1e1211SMatthew G Knepley   */
252cb1e1211SMatthew G Knepley 
25376185916SToby Isaac   ierr = DMPlexComputeAnchorAdjacencies(dm,section,sectionGlobal,&anchorSectionAdj,&anchorAdj);CHKERRQ(ierr);
25476185916SToby Isaac 
255cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
25676185916SToby Isaac     PetscInt dof, off, d, q, anDof;
25770034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
258cb1e1211SMatthew G Knepley 
259fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
260cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
261cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
26270034214SMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr);
263cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
264cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
265cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
266cb1e1211SMatthew G Knepley 
267cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
268cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
269cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
270cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
271cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
272cb1e1211SMatthew G Knepley       }
273cb1e1211SMatthew G Knepley     }
27476185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
27576185916SToby Isaac     if (anDof) {
27676185916SToby Isaac       for (d = off; d < off+dof; ++d) {
27776185916SToby Isaac         ierr = PetscSectionAddDof(leafSectionAdj, d, anDof);CHKERRQ(ierr);
27876185916SToby Isaac       }
27976185916SToby Isaac     }
280cb1e1211SMatthew G Knepley   }
281cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
282cb1e1211SMatthew G Knepley   if (debug) {
283cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n");CHKERRQ(ierr);
284cb1e1211SMatthew G Knepley     ierr = PetscSectionView(leafSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
285cb1e1211SMatthew G Knepley   }
286cb1e1211SMatthew G Knepley   /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
28747a6104aSMatthew G. Knepley   if (doComm) {
288cb1e1211SMatthew G Knepley     ierr = PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
289cb1e1211SMatthew G Knepley     ierr = PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
290cb1e1211SMatthew G Knepley   }
291cb1e1211SMatthew G Knepley   if (debug) {
292cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n");CHKERRQ(ierr);
293cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
294cb1e1211SMatthew G Knepley   }
295cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
296cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
29776185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
298cb1e1211SMatthew G Knepley 
299cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
300cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
301cb1e1211SMatthew G Knepley     if (!dof) continue;
302cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
303cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
30470034214SMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr);
305cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
306cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
307cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
308cb1e1211SMatthew G Knepley 
309cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
310cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
311cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
312cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
313cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
314cb1e1211SMatthew G Knepley       }
315cb1e1211SMatthew G Knepley     }
31676185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
31776185916SToby Isaac     if (anDof) {
31876185916SToby Isaac       for (d = off; d < off+dof; ++d) {
31976185916SToby Isaac         ierr = PetscSectionAddDof(rootSectionAdj, d, anDof);CHKERRQ(ierr);
32076185916SToby Isaac       }
32176185916SToby Isaac     }
322cb1e1211SMatthew G Knepley   }
323cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
324cb1e1211SMatthew G Knepley   if (debug) {
325cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n");CHKERRQ(ierr);
326cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
327cb1e1211SMatthew G Knepley   }
328cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
329cb1e1211SMatthew G Knepley   ierr = PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets);CHKERRQ(ierr);
330cb1e1211SMatthew G Knepley   ierr = PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj);CHKERRQ(ierr);
331cb1e1211SMatthew G Knepley   if (debug) {
332cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency SF for Preallocation:\n");CHKERRQ(ierr);
333cb1e1211SMatthew G Knepley     ierr = PetscSFView(sfAdj, NULL);CHKERRQ(ierr);
334cb1e1211SMatthew G Knepley   }
335cb1e1211SMatthew G Knepley   ierr = PetscSFDestroy(&sfDof);CHKERRQ(ierr);
336cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
337cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
338cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(leafSectionAdj, &adjSize);CHKERRQ(ierr);
3391795a4d1SJed Brown   ierr = PetscCalloc1(adjSize, &adj);CHKERRQ(ierr);
340cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
34176185916SToby Isaac     PetscInt dof, off, d, q, anDof, anOff;
34270034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
343cb1e1211SMatthew G Knepley 
344fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
345cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
346cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
34770034214SMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr);
34876185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
34976185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
350cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
351cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
352cb1e1211SMatthew G Knepley 
353cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(leafSectionAdj, d, &aoff);CHKERRQ(ierr);
354cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
355cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
356cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
357cb1e1211SMatthew G Knepley 
358cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
359cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
360cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
361cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
362cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
363cb1e1211SMatthew G Knepley           adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd;
364cb1e1211SMatthew G Knepley           ++i;
365cb1e1211SMatthew G Knepley         }
366cb1e1211SMatthew G Knepley       }
36776185916SToby Isaac       for (q = 0; q < anDof; q++) {
36876185916SToby Isaac         adj[aoff+i] = anchorAdj[anOff+q];
36976185916SToby Isaac         ++i;
37076185916SToby Isaac       }
371cb1e1211SMatthew G Knepley     }
372cb1e1211SMatthew G Knepley   }
373cb1e1211SMatthew G Knepley   /* Debugging */
374cb1e1211SMatthew G Knepley   if (debug) {
375cb1e1211SMatthew G Knepley     IS tmp;
376cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Leaf adjacency indices\n");CHKERRQ(ierr);
377cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
378cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
3790a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
380cb1e1211SMatthew G Knepley   }
381*543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
382cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(rootSectionAdj, &adjSize);CHKERRQ(ierr);
383785e854fSJed Brown   ierr = PetscMalloc1(adjSize, &rootAdj);CHKERRQ(ierr);
384cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
38547a6104aSMatthew G. Knepley   if (doComm) {
386*543482b8SMatthew G. Knepley     const PetscInt *indegree;
387*543482b8SMatthew G. Knepley     PetscInt       *remoteadj, radjsize = 0;
388*543482b8SMatthew G. Knepley 
389*543482b8SMatthew G. Knepley     ierr = PetscSFComputeDegreeBegin(sfAdj, &indegree);CHKERRQ(ierr);
390*543482b8SMatthew G. Knepley     ierr = PetscSFComputeDegreeEnd(sfAdj, &indegree);CHKERRQ(ierr);
391*543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
392*543482b8SMatthew G. Knepley     ierr = PetscMalloc1(radjsize, &remoteadj);CHKERRQ(ierr);
393*543482b8SMatthew G. Knepley     ierr = PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr);
394*543482b8SMatthew G. Knepley     ierr = PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr);
395*543482b8SMatthew G. Knepley     for (p = 0, r = 0; p < adjSize; ++p) {
396*543482b8SMatthew G. Knepley       PetscInt s;
397*543482b8SMatthew G. Knepley       for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[p] = remoteadj[r];
398*543482b8SMatthew G. Knepley     }
399*543482b8SMatthew G. Knepley     if (r != radjsize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", r, radjsize);
400*543482b8SMatthew G. Knepley     ierr = PetscFree(remoteadj);CHKERRQ(ierr);
401cb1e1211SMatthew G Knepley   }
402cb1e1211SMatthew G Knepley   ierr = PetscSFDestroy(&sfAdj);CHKERRQ(ierr);
403cb1e1211SMatthew G Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
404cb1e1211SMatthew G Knepley   /* Debugging */
405cb1e1211SMatthew G Knepley   if (debug) {
406cb1e1211SMatthew G Knepley     IS tmp;
407cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after gather\n");CHKERRQ(ierr);
408cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
409cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
4100a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
411cb1e1211SMatthew G Knepley   }
412cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
413cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
41476185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
415cb1e1211SMatthew G Knepley 
416cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
417cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
418cb1e1211SMatthew G Knepley     if (!dof) continue;
419cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
420cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
42170034214SMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr);
42276185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
42376185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
424cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
425cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
426cb1e1211SMatthew G Knepley 
427cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
428cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
429cb1e1211SMatthew G Knepley       i    = adof-1;
43076185916SToby Isaac       for (q = 0; q < anDof; q++) {
43176185916SToby Isaac         rootAdj[aoff+i] = anchorAdj[anOff+q];
43276185916SToby Isaac         --i;
43376185916SToby Isaac       }
434cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
435cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
436cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
437cb1e1211SMatthew G Knepley 
438cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
439cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
440cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
441cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
442cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
443cb1e1211SMatthew G Knepley           rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
444cb1e1211SMatthew G Knepley           --i;
445cb1e1211SMatthew G Knepley         }
446cb1e1211SMatthew G Knepley       }
447cb1e1211SMatthew G Knepley     }
448cb1e1211SMatthew G Knepley   }
449cb1e1211SMatthew G Knepley   /* Debugging */
450cb1e1211SMatthew G Knepley   if (debug) {
451cb1e1211SMatthew G Knepley     IS tmp;
452cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices\n");CHKERRQ(ierr);
453cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
454cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
4550a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
456cb1e1211SMatthew G Knepley   }
457cb1e1211SMatthew G Knepley   /* Compress indices */
458cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
459cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
460cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
461cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
462cb1e1211SMatthew G Knepley 
463cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
464cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
465cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
466cb1e1211SMatthew G Knepley     if (!dof) continue;
467cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
468cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
469cb1e1211SMatthew G Knepley     for (d = off; d < off+dof-cdof; ++d) {
470cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
471cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
472cb1e1211SMatthew G Knepley       ierr = PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]);CHKERRQ(ierr);
473cb1e1211SMatthew G Knepley       ierr = PetscSectionSetDof(rootSectionAdj, d, adof);CHKERRQ(ierr);
474cb1e1211SMatthew G Knepley     }
475cb1e1211SMatthew G Knepley   }
476cb1e1211SMatthew G Knepley   /* Debugging */
477cb1e1211SMatthew G Knepley   if (debug) {
478cb1e1211SMatthew G Knepley     IS tmp;
479cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n");CHKERRQ(ierr);
480cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
481cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after compression\n");CHKERRQ(ierr);
482cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
483cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
4840a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
485cb1e1211SMatthew G Knepley   }
486cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
487cb1e1211SMatthew G Knepley   ierr = PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd);CHKERRQ(ierr);
488cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &sectionAdj);CHKERRQ(ierr);
489cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd);CHKERRQ(ierr);
490cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
49176185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
492cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
493cb1e1211SMatthew G Knepley 
494cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
495cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
496cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
497cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
498cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
499cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
500cb1e1211SMatthew G Knepley 
501cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
502cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
503cb1e1211SMatthew G Knepley       if (ldof > 0) {
504cb1e1211SMatthew G Knepley         /* We do not own this point */
505cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
506cb1e1211SMatthew G Knepley         ierr = PetscSectionSetDof(sectionAdj, goff+d, rdof);CHKERRQ(ierr);
507cb1e1211SMatthew G Knepley       } else {
508cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
509cb1e1211SMatthew G Knepley       }
510cb1e1211SMatthew G Knepley     }
511cb1e1211SMatthew G Knepley     if (found) continue;
512cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
513cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
51470034214SMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr);
515cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
516cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
517cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof, noff;
518cb1e1211SMatthew G Knepley 
519cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
520cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
521cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
522cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(section, padj, &noff);CHKERRQ(ierr);
523cb1e1211SMatthew G Knepley       for (d = goff; d < goff+dof-cdof; ++d) {
524cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(sectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
525cb1e1211SMatthew G Knepley       }
526cb1e1211SMatthew G Knepley     }
52776185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
52876185916SToby Isaac     if (anDof) {
52976185916SToby Isaac       for (d = goff; d < goff+dof-cdof; ++d) {
53076185916SToby Isaac         ierr = PetscSectionAddDof(sectionAdj, d, anDof);CHKERRQ(ierr);
53176185916SToby Isaac       }
53276185916SToby Isaac     }
533cb1e1211SMatthew G Knepley   }
534cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(sectionAdj);CHKERRQ(ierr);
535cb1e1211SMatthew G Knepley   if (debug) {
536cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation:\n");CHKERRQ(ierr);
537cb1e1211SMatthew G Knepley     ierr = PetscSectionView(sectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
538cb1e1211SMatthew G Knepley   }
539cb1e1211SMatthew G Knepley   /* Get adjacent indices */
540cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(sectionAdj, &numCols);CHKERRQ(ierr);
541785e854fSJed Brown   ierr = PetscMalloc1(numCols, &cols);CHKERRQ(ierr);
542cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
54376185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
544cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
545cb1e1211SMatthew G Knepley 
546cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
547cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
548cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
549cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
550cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
551cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
552cb1e1211SMatthew G Knepley 
553cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
554cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
555cb1e1211SMatthew G Knepley       if (ldof > 0) {
556cb1e1211SMatthew G Knepley         /* We do not own this point */
557cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
558cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
559cb1e1211SMatthew G Knepley 
560cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionAdj, goff+d, &aoff);CHKERRQ(ierr);
561cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(rootSectionAdj, off+d, &roff);CHKERRQ(ierr);
562cb1e1211SMatthew G Knepley         ierr = PetscMemcpy(&cols[aoff], &rootAdj[roff], rdof * sizeof(PetscInt));CHKERRQ(ierr);
563cb1e1211SMatthew G Knepley       } else {
564cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
565cb1e1211SMatthew G Knepley       }
566cb1e1211SMatthew G Knepley     }
567cb1e1211SMatthew G Knepley     if (found) continue;
56870034214SMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr);
56976185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
57076185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
571cb1e1211SMatthew G Knepley     for (d = goff; d < goff+dof-cdof; ++d) {
572cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
573cb1e1211SMatthew G Knepley 
574cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(sectionAdj, d, &adof);CHKERRQ(ierr);
575cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(sectionAdj, d, &aoff);CHKERRQ(ierr);
576cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
577cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
578cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
579cb1e1211SMatthew G Knepley         const PetscInt *ncind;
580cb1e1211SMatthew G Knepley 
581cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
582cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
583cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
584cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
585cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintIndices(section, padj, &ncind);CHKERRQ(ierr);
586cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
587cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd, ++i) {
588cb1e1211SMatthew G Knepley           cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
589cb1e1211SMatthew G Knepley         }
590cb1e1211SMatthew G Knepley       }
59176185916SToby Isaac       for (q = 0; q < anDof; q++, i++) {
59276185916SToby Isaac         cols[aoff+i] = anchorAdj[anOff + q];
59376185916SToby Isaac       }
594cb1e1211SMatthew G Knepley       if (i != adof) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of entries %D != %D for dof %D (point %D)", i, adof, d, p);
595cb1e1211SMatthew G Knepley     }
596cb1e1211SMatthew G Knepley   }
59776185916SToby Isaac   ierr = PetscSectionDestroy(&anchorSectionAdj);CHKERRQ(ierr);
598cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&leafSectionAdj);CHKERRQ(ierr);
599cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&rootSectionAdj);CHKERRQ(ierr);
60076185916SToby Isaac   ierr = PetscFree(anchorAdj);CHKERRQ(ierr);
601cb1e1211SMatthew G Knepley   ierr = PetscFree(rootAdj);CHKERRQ(ierr);
60270034214SMatthew G. Knepley   ierr = PetscFree(tmpAdj);CHKERRQ(ierr);
603cb1e1211SMatthew G Knepley   /* Debugging */
604cb1e1211SMatthew G Knepley   if (debug) {
605cb1e1211SMatthew G Knepley     IS tmp;
606cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Column indices\n");CHKERRQ(ierr);
607cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
608cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
6090a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
610cb1e1211SMatthew G Knepley   }
611cb1e1211SMatthew G Knepley   /* Create allocation vectors from adjacency graph */
612cb1e1211SMatthew G Knepley   ierr = MatGetLocalSize(A, &locRows, NULL);CHKERRQ(ierr);
613cb1e1211SMatthew G Knepley   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)A), &rLayout);CHKERRQ(ierr);
614cb1e1211SMatthew G Knepley   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
615cb1e1211SMatthew G Knepley   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
616cb1e1211SMatthew G Knepley   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
617cb1e1211SMatthew G Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
618cb1e1211SMatthew G Knepley   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
619cb1e1211SMatthew G Knepley   /* Only loop over blocks of rows */
620cb1e1211SMatthew G Knepley   if (rStart%bs || rEnd%bs) SETERRQ3(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Invalid layout [%d, %d) for matrix, must be divisible by block size %d", rStart, rEnd, bs);
621cb1e1211SMatthew G Knepley   for (r = rStart/bs; r < rEnd/bs; ++r) {
622cb1e1211SMatthew G Knepley     const PetscInt row = r*bs;
623cb1e1211SMatthew G Knepley     PetscInt       numCols, cStart, c;
624cb1e1211SMatthew G Knepley 
625cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(sectionAdj, row, &numCols);CHKERRQ(ierr);
626cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionAdj, row, &cStart);CHKERRQ(ierr);
627cb1e1211SMatthew G Knepley     for (c = cStart; c < cStart+numCols; ++c) {
628cb1e1211SMatthew G Knepley       if ((cols[c] >= rStart*bs) && (cols[c] < rEnd*bs)) {
629cb1e1211SMatthew G Knepley         ++dnz[r-rStart];
630cb1e1211SMatthew G Knepley         if (cols[c] >= row) ++dnzu[r-rStart];
631cb1e1211SMatthew G Knepley       } else {
632cb1e1211SMatthew G Knepley         ++onz[r-rStart];
633cb1e1211SMatthew G Knepley         if (cols[c] >= row) ++onzu[r-rStart];
634cb1e1211SMatthew G Knepley       }
635cb1e1211SMatthew G Knepley     }
636cb1e1211SMatthew G Knepley   }
637cb1e1211SMatthew G Knepley   if (bs > 1) {
638cb1e1211SMatthew G Knepley     for (r = 0; r < locRows/bs; ++r) {
639cb1e1211SMatthew G Knepley       dnz[r]  /= bs;
640cb1e1211SMatthew G Knepley       onz[r]  /= bs;
641cb1e1211SMatthew G Knepley       dnzu[r] /= bs;
642cb1e1211SMatthew G Knepley       onzu[r] /= bs;
643cb1e1211SMatthew G Knepley     }
644cb1e1211SMatthew G Knepley   }
645cb1e1211SMatthew G Knepley   /* Set matrix pattern */
646cb1e1211SMatthew G Knepley   ierr = MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu);CHKERRQ(ierr);
647cb1e1211SMatthew G Knepley   ierr = MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
64889545effSMatthew G. Knepley   /* Check for symmetric storage */
64989545effSMatthew G. Knepley   ierr = MatGetType(A, &mtype);CHKERRQ(ierr);
65089545effSMatthew G. Knepley   ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr);
65189545effSMatthew G. Knepley   ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr);
65289545effSMatthew G. Knepley   ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr);
65389545effSMatthew G. Knepley   if (isSymBlock || isSymSeqBlock || isSymMPIBlock) {ierr = MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);}
654cb1e1211SMatthew G Knepley   /* Fill matrix with zeros */
655cb1e1211SMatthew G Knepley   if (fillMatrix) {
656cb1e1211SMatthew G Knepley     PetscScalar *values;
657cb1e1211SMatthew G Knepley     PetscInt     maxRowLen = 0;
658cb1e1211SMatthew G Knepley 
659cb1e1211SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
660cb1e1211SMatthew G Knepley       PetscInt len;
661cb1e1211SMatthew G Knepley 
662cb1e1211SMatthew G Knepley       ierr      = PetscSectionGetDof(sectionAdj, r, &len);CHKERRQ(ierr);
663cb1e1211SMatthew G Knepley       maxRowLen = PetscMax(maxRowLen, len);
664cb1e1211SMatthew G Knepley     }
6651795a4d1SJed Brown     ierr = PetscCalloc1(maxRowLen, &values);CHKERRQ(ierr);
666cb1e1211SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
667cb1e1211SMatthew G Knepley       PetscInt numCols, cStart;
668cb1e1211SMatthew G Knepley 
669cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
670cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
671cb1e1211SMatthew G Knepley       ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr);
672cb1e1211SMatthew G Knepley     }
673cb1e1211SMatthew G Knepley     ierr = PetscFree(values);CHKERRQ(ierr);
674cb1e1211SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
675cb1e1211SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
676cb1e1211SMatthew G Knepley   }
677e228b242SToby Isaac   /* restore original useAnchors */
678e228b242SToby Isaac   ierr = DMPlexSetAdjacencyUseAnchors(dm,useAnchors);CHKERRQ(ierr);
679cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&sectionAdj);CHKERRQ(ierr);
680cb1e1211SMatthew G Knepley   ierr = PetscFree(cols);CHKERRQ(ierr);
681a72f3261SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr);
682cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
683cb1e1211SMatthew G Knepley }
684cb1e1211SMatthew G Knepley 
685cb1e1211SMatthew G Knepley #if 0
686cb1e1211SMatthew G Knepley #undef __FUNCT__
687cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator_2"
688cb1e1211SMatthew 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)
689cb1e1211SMatthew G Knepley {
690cb1e1211SMatthew G Knepley   PetscInt       *tmpClosure,*tmpAdj,*visits;
691cb1e1211SMatthew G Knepley   PetscInt        c,cStart,cEnd,pStart,pEnd;
692cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
693cb1e1211SMatthew G Knepley 
694cb1e1211SMatthew G Knepley   PetscFunctionBegin;
695c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
696cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
697cb1e1211SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
698cb1e1211SMatthew G Knepley 
699e7b80042SMatthew G. Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1));
700cb1e1211SMatthew G Knepley 
701cb1e1211SMatthew G Knepley   ierr    = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
702cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
703cb1e1211SMatthew G Knepley 
704dcca6d9dSJed Brown   ierr = PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits);CHKERRQ(ierr);
705cb1e1211SMatthew G Knepley   ierr = PetscMemzero(lvisits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr);
706cb1e1211SMatthew G Knepley   ierr = PetscMemzero(visits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr);
707cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
708cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
709cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
710cb1e1211SMatthew G Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support);CHKERRQ(ierr);
711cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
712cb1e1211SMatthew G Knepley   }
713cb1e1211SMatthew G Knepley   ierr = PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr);
714cb1e1211SMatthew G Knepley   ierr = PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr);
715cb1e1211SMatthew G Knepley   ierr = PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr);
716cb1e1211SMatthew G Knepley   ierr = PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr);
717cb1e1211SMatthew G Knepley 
718cb1e1211SMatthew G Knepley   ierr = PetscSFGetRanks();CHKERRQ(ierr);
719cb1e1211SMatthew G Knepley 
720cb1e1211SMatthew G Knepley 
721dcca6d9dSJed Brown   ierr = PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner);CHKERRQ(ierr);
722cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
723cb1e1211SMatthew G Knepley     ierr = PetscMemzero(cellmat,maxClosureSize*maxClosureSize*sizeof(PetscInt));CHKERRQ(ierr);
724cb1e1211SMatthew G Knepley     /*
725cb1e1211SMatthew G Knepley      Depth-first walk of transitive closure.
726cb1e1211SMatthew 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.
727cb1e1211SMatthew G Knepley      This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
728cb1e1211SMatthew G Knepley      */
729cb1e1211SMatthew G Knepley   }
730cb1e1211SMatthew G Knepley 
731cb1e1211SMatthew G Knepley   ierr = PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM);CHKERRQ(ierr);
732cb1e1211SMatthew G Knepley   ierr = PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM);CHKERRQ(ierr);
733cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
734cb1e1211SMatthew G Knepley }
735cb1e1211SMatthew G Knepley #endif
736