xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision e7bcfa36292a613e38c1fdf696c3f2f4aa82e031)
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 
6cb1e1211SMatthew G Knepley #undef __FUNCT__
776185916SToby Isaac #define __FUNCT__ "DMPlexComputeAnchorAdjacencies"
876185916SToby Isaac /* get adjacencies due to point-to-point constraints that can't be found with DMPlexGetAdjacency() */
9e101f074SMatthew G. Knepley static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[])
1076185916SToby Isaac {
1176185916SToby Isaac   PetscInt       pStart, pEnd;
12e101f074SMatthew G. Knepley   PetscSection   section, sectionGlobal, adjSec, aSec;
1376185916SToby Isaac   IS             aIS;
1476185916SToby Isaac   PetscErrorCode ierr;
1576185916SToby Isaac 
1676185916SToby Isaac   PetscFunctionBegin;
17e101f074SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
18e101f074SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
1976185916SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) section), &adjSec);CHKERRQ(ierr);
2076185916SToby Isaac   ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr);
2176185916SToby Isaac   ierr = PetscSectionSetChart(adjSec,pStart,pEnd);CHKERRQ(ierr);
2276185916SToby Isaac 
23a17985deSToby Isaac   ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
2476185916SToby Isaac   if (aSec) {
2576185916SToby Isaac     const PetscInt *anchors;
2676185916SToby Isaac     PetscInt       p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize;
2776185916SToby Isaac     PetscInt       *tmpAdjP = NULL, *tmpAdjQ = NULL;
2876185916SToby Isaac     PetscSection   inverseSec;
2976185916SToby Isaac 
3076185916SToby Isaac     /* invert the constraint-to-anchor map */
3176185916SToby Isaac     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)aSec),&inverseSec);CHKERRQ(ierr);
3276185916SToby Isaac     ierr = PetscSectionSetChart(inverseSec,pStart,pEnd);CHKERRQ(ierr);
3376185916SToby Isaac     ierr = ISGetLocalSize(aIS, &aSize);CHKERRQ(ierr);
3476185916SToby Isaac     ierr = ISGetIndices(aIS, &anchors);CHKERRQ(ierr);
3576185916SToby Isaac 
3676185916SToby Isaac     for (p = 0; p < aSize; p++) {
3776185916SToby Isaac       PetscInt a = anchors[p];
3876185916SToby Isaac 
3976185916SToby Isaac       ierr = PetscSectionAddDof(inverseSec,a,1);CHKERRQ(ierr);
4076185916SToby Isaac     }
4176185916SToby Isaac     ierr = PetscSectionSetUp(inverseSec);CHKERRQ(ierr);
4276185916SToby Isaac     ierr = PetscSectionGetStorageSize(inverseSec,&iSize);CHKERRQ(ierr);
4376185916SToby Isaac     ierr = PetscMalloc1(iSize,&inverse);CHKERRQ(ierr);
4476185916SToby Isaac     ierr = PetscCalloc1(pEnd-pStart,&offsets);CHKERRQ(ierr);
4576185916SToby Isaac     ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
4676185916SToby Isaac     for (p = aStart; p < aEnd; p++) {
4776185916SToby Isaac       PetscInt dof, off;
4876185916SToby Isaac 
4976185916SToby Isaac       ierr = PetscSectionGetDof(aSec, p, &dof);CHKERRQ(ierr);
5076185916SToby Isaac       ierr = PetscSectionGetOffset(aSec, p, &off);CHKERRQ(ierr);
5176185916SToby Isaac 
5276185916SToby Isaac       for (q = 0; q < dof; q++) {
5376185916SToby Isaac         PetscInt iOff;
5476185916SToby Isaac 
5576185916SToby Isaac         a = anchors[off + q];
5676185916SToby Isaac         ierr = PetscSectionGetOffset(inverseSec, a, &iOff);CHKERRQ(ierr);
5776185916SToby Isaac         inverse[iOff + offsets[a-pStart]++] = p;
5876185916SToby Isaac       }
5976185916SToby Isaac     }
6076185916SToby Isaac     ierr = ISRestoreIndices(aIS, &anchors);CHKERRQ(ierr);
6176185916SToby Isaac     ierr = PetscFree(offsets);CHKERRQ(ierr);
6276185916SToby Isaac 
6376185916SToby Isaac     /* construct anchorAdj and adjSec
6476185916SToby Isaac      *
6576185916SToby Isaac      * loop over anchors:
6676185916SToby Isaac      *   construct anchor adjacency
6776185916SToby Isaac      *   loop over constrained:
6876185916SToby Isaac      *     construct constrained adjacency
6976185916SToby Isaac      *     if not in anchor adjacency, add to dofs
7076185916SToby Isaac      * setup adjSec, allocate anchorAdj
7176185916SToby Isaac      * loop over anchors:
7276185916SToby Isaac      *   construct anchor adjacency
7376185916SToby Isaac      *   loop over constrained:
7476185916SToby Isaac      *     construct constrained adjacency
7576185916SToby Isaac      *     if not in anchor adjacency
7676185916SToby Isaac      *       if not already in list, put in list
7776185916SToby Isaac      *   sort, unique, reduce dof count
7876185916SToby Isaac      * optional: compactify
7976185916SToby Isaac      */
8076185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
8176185916SToby Isaac       PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE;
8276185916SToby Isaac 
8376185916SToby Isaac       ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr);
8476185916SToby Isaac       if (!iDof) continue;
8576185916SToby Isaac       ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr);
86a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP);CHKERRQ(ierr);
8776185916SToby Isaac       for (i = 0; i < iDof; i++) {
8876185916SToby Isaac         PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE;
8976185916SToby Isaac 
9076185916SToby Isaac         q = inverse[iOff + i];
91a9fb6443SMatthew G. Knepley         ierr = DMPlexGetAdjacency_Internal(dm,q,useCone,useClosure,PETSC_TRUE,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr);
9276185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
9376185916SToby Isaac           qAdj = tmpAdjQ[r];
9476185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
9576185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
9676185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
9776185916SToby Isaac           }
9876185916SToby Isaac           if (s < numAdjP) continue;
9976185916SToby Isaac           ierr  = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr);
10076185916SToby Isaac           ierr  = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr);
10176185916SToby Isaac           iNew += qAdjDof - qAdjCDof;
10276185916SToby Isaac         }
10376185916SToby Isaac         ierr = PetscSectionAddDof(adjSec,p,iNew);CHKERRQ(ierr);
10476185916SToby Isaac       }
10576185916SToby Isaac     }
10676185916SToby Isaac 
10776185916SToby Isaac     ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr);
10876185916SToby Isaac     ierr = PetscSectionGetStorageSize(adjSec,&adjSize);CHKERRQ(ierr);
10976185916SToby Isaac     ierr = PetscMalloc1(adjSize,&adj);CHKERRQ(ierr);
11076185916SToby Isaac 
11176185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
11276185916SToby Isaac       PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE;
11376185916SToby Isaac 
11476185916SToby Isaac       ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr);
11576185916SToby Isaac       if (!iDof) continue;
11676185916SToby Isaac       ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr);
117a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP);CHKERRQ(ierr);
11876185916SToby Isaac       ierr = PetscSectionGetDof(adjSec,p,&aDof);CHKERRQ(ierr);
11976185916SToby Isaac       ierr = PetscSectionGetOffset(adjSec,p,&aOff);CHKERRQ(ierr);
12076185916SToby Isaac       aOffOrig = aOff;
12176185916SToby Isaac       for (i = 0; i < iDof; i++) {
12276185916SToby Isaac         PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE;
12376185916SToby Isaac 
12476185916SToby Isaac         q = inverse[iOff + i];
125a9fb6443SMatthew G. Knepley         ierr = DMPlexGetAdjacency_Internal(dm,q,useCone,useClosure,PETSC_TRUE,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr);
12676185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
12776185916SToby Isaac           qAdj = tmpAdjQ[r];
12876185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
12976185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
13076185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
13176185916SToby Isaac           }
13276185916SToby Isaac           if (s < numAdjP) continue;
13376185916SToby Isaac           ierr  = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr);
13476185916SToby Isaac           ierr  = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr);
13576185916SToby Isaac           ierr  = PetscSectionGetOffset(sectionGlobal,qAdj,&qAdjOff);CHKERRQ(ierr);
13676185916SToby Isaac           for (nd = 0; nd < qAdjDof-qAdjCDof; ++nd) {
13776185916SToby Isaac             adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff+1) : qAdjOff) + nd;
13876185916SToby Isaac           }
13976185916SToby Isaac         }
14076185916SToby Isaac       }
14176185916SToby Isaac       ierr = PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]);CHKERRQ(ierr);
14276185916SToby Isaac       ierr = PetscSectionSetDof(adjSec,p,aDof);CHKERRQ(ierr);
14376185916SToby Isaac     }
14476185916SToby Isaac     *anchorAdj = adj;
14576185916SToby Isaac 
14676185916SToby Isaac     /* clean up */
14776185916SToby Isaac     ierr = PetscSectionDestroy(&inverseSec);CHKERRQ(ierr);
14876185916SToby Isaac     ierr = PetscFree(inverse);CHKERRQ(ierr);
14976185916SToby Isaac     ierr = PetscFree(tmpAdjP);CHKERRQ(ierr);
15076185916SToby Isaac     ierr = PetscFree(tmpAdjQ);CHKERRQ(ierr);
15176185916SToby Isaac   }
15276185916SToby Isaac   else {
15376185916SToby Isaac     *anchorAdj = NULL;
15476185916SToby Isaac     ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr);
15576185916SToby Isaac   }
15676185916SToby Isaac   *anchorSectionAdj = adjSec;
15776185916SToby Isaac   PetscFunctionReturn(0);
15876185916SToby Isaac }
15976185916SToby Isaac 
16076185916SToby Isaac #undef __FUNCT__
161a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateAdjacencySection_Static"
162e101f074SMatthew G. Knepley PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx)
163cb1e1211SMatthew G Knepley {
164cb1e1211SMatthew G Knepley   MPI_Comm           comm;
165a9fb6443SMatthew G. Knepley   PetscMPIInt        size;
166a9fb6443SMatthew G. Knepley   PetscBool          doCommLocal, doComm, debug = PETSC_FALSE;
167a9fb6443SMatthew G. Knepley   PetscSF            sf, sfAdj;
168e101f074SMatthew G. Knepley   PetscSection       section, sectionGlobal, leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj;
169a9fb6443SMatthew G. Knepley   PetscInt           nroots, nleaves, l, p, r;
170cb1e1211SMatthew G Knepley   const PetscInt    *leaves;
171cb1e1211SMatthew G Knepley   const PetscSFNode *remotes;
172cb1e1211SMatthew G Knepley   PetscInt           dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols;
1738821704fSMatthew G. Knepley   PetscInt          *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets;
17470034214SMatthew G. Knepley   PetscInt           adjSize;
175cb1e1211SMatthew G Knepley   PetscErrorCode     ierr;
176cb1e1211SMatthew G Knepley 
177cb1e1211SMatthew G Knepley   PetscFunctionBegin;
178cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
179c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr);
180cb1e1211SMatthew G Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
181c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
182cb1e1211SMatthew G Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
183e101f074SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
184e101f074SMatthew G. Knepley   ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
18547a6104aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
18647a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
187b2566f29SBarry Smith   ierr = MPIU_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr);
188cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
189cb1e1211SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
190cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &numDof);CHKERRQ(ierr);
191cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &leafSectionAdj);CHKERRQ(ierr);
192cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(leafSectionAdj, 0, numDof);CHKERRQ(ierr);
193cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &rootSectionAdj);CHKERRQ(ierr);
194cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(rootSectionAdj, 0, numDof);CHKERRQ(ierr);
195cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
196cb1e1211SMatthew G Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes);CHKERRQ(ierr);
197cb1e1211SMatthew G Knepley   /*
198964bf7afSMatthew G. Knepley    section        - maps points to (# dofs, local dofs)
199964bf7afSMatthew G. Knepley    sectionGlobal  - maps points to (# dofs, global dofs)
200964bf7afSMatthew G. Knepley    leafSectionAdj - maps unowned local dofs to # adj dofs
201964bf7afSMatthew G. Knepley    rootSectionAdj - maps   owned local dofs to # adj dofs
202964bf7afSMatthew G. Knepley    adj            - adj global dofs indexed by leafSectionAdj
203964bf7afSMatthew G. Knepley    rootAdj        - adj global dofs indexed by rootSectionAdj
204964bf7afSMatthew G. Knepley    sf    - describes shared points across procs
205964bf7afSMatthew G. Knepley    sfDof - describes shared dofs across procs
206964bf7afSMatthew G. Knepley    sfAdj - describes shared adjacent dofs across procs
207cb1e1211SMatthew G Knepley    ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
20876185916SToby Isaac   (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj
20976185916SToby Isaac        (This is done in DMPlexComputeAnchorAdjacencies())
210cb1e1211SMatthew G Knepley     1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
211cb1e1211SMatthew G Knepley        Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
212cb1e1211SMatthew G Knepley     2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
213cb1e1211SMatthew G Knepley        Create sfAdj connecting rootSectionAdj and leafSectionAdj
214cb1e1211SMatthew G Knepley     3. Visit unowned points on interface, write adjacencies to adj
215cb1e1211SMatthew G Knepley        Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
216cb1e1211SMatthew G Knepley     4. Visit owned points on interface, write adjacencies to rootAdj
217cb1e1211SMatthew G Knepley        Remove redundancy in rootAdj
218cb1e1211SMatthew G Knepley    ** The last two traversals use transitive closure
219cb1e1211SMatthew G Knepley     5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
220cb1e1211SMatthew G Knepley        Allocate memory addressed by sectionAdj (cols)
221cb1e1211SMatthew G Knepley     6. Visit all owned points in the subdomain, insert dof adjacencies into cols
222cb1e1211SMatthew G Knepley    ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
223cb1e1211SMatthew G Knepley   */
224e101f074SMatthew G. Knepley   ierr = DMPlexComputeAnchorAdjacencies(dm, useCone, useClosure, &anchorSectionAdj, &anchorAdj);CHKERRQ(ierr);
225cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
22676185916SToby Isaac     PetscInt dof, off, d, q, anDof;
22770034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
228cb1e1211SMatthew G Knepley 
229fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
230cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
231cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
232a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
233cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
234cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
235cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
236cb1e1211SMatthew G Knepley 
237cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
238cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
239cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
240cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
241cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
242cb1e1211SMatthew G Knepley       }
243cb1e1211SMatthew G Knepley     }
24476185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
24576185916SToby Isaac     if (anDof) {
24676185916SToby Isaac       for (d = off; d < off+dof; ++d) {
24776185916SToby Isaac         ierr = PetscSectionAddDof(leafSectionAdj, d, anDof);CHKERRQ(ierr);
24876185916SToby Isaac       }
24976185916SToby Isaac     }
250cb1e1211SMatthew G Knepley   }
251cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
252cb1e1211SMatthew G Knepley   if (debug) {
253cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n");CHKERRQ(ierr);
254cb1e1211SMatthew G Knepley     ierr = PetscSectionView(leafSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
255cb1e1211SMatthew G Knepley   }
256cb1e1211SMatthew G Knepley   /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
25747a6104aSMatthew G. Knepley   if (doComm) {
258cb1e1211SMatthew G Knepley     ierr = PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
259cb1e1211SMatthew G Knepley     ierr = PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
260cb1e1211SMatthew G Knepley   }
261cb1e1211SMatthew G Knepley   if (debug) {
262cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n");CHKERRQ(ierr);
263cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
264cb1e1211SMatthew G Knepley   }
265cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
266cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
26776185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
268cb1e1211SMatthew G Knepley 
269cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
270cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
271cb1e1211SMatthew G Knepley     if (!dof) continue;
272cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
273cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
274a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
275cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
276cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
277cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
278cb1e1211SMatthew G Knepley 
279cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
280cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
281cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
282cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
283cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
284cb1e1211SMatthew G Knepley       }
285cb1e1211SMatthew G Knepley     }
28676185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
28776185916SToby Isaac     if (anDof) {
28876185916SToby Isaac       for (d = off; d < off+dof; ++d) {
28976185916SToby Isaac         ierr = PetscSectionAddDof(rootSectionAdj, d, anDof);CHKERRQ(ierr);
29076185916SToby Isaac       }
29176185916SToby Isaac     }
292cb1e1211SMatthew G Knepley   }
293cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
294cb1e1211SMatthew G Knepley   if (debug) {
295cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n");CHKERRQ(ierr);
296cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
297cb1e1211SMatthew G Knepley   }
298cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
299cb1e1211SMatthew G Knepley   ierr = PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets);CHKERRQ(ierr);
300cb1e1211SMatthew G Knepley   ierr = PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj);CHKERRQ(ierr);
3018a713f3bSLawrence Mitchell   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
302cb1e1211SMatthew G Knepley   if (debug) {
303cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency SF for Preallocation:\n");CHKERRQ(ierr);
304cb1e1211SMatthew G Knepley     ierr = PetscSFView(sfAdj, NULL);CHKERRQ(ierr);
305cb1e1211SMatthew G Knepley   }
306cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
307cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
308cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(leafSectionAdj, &adjSize);CHKERRQ(ierr);
3091795a4d1SJed Brown   ierr = PetscCalloc1(adjSize, &adj);CHKERRQ(ierr);
310cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
31176185916SToby Isaac     PetscInt dof, off, d, q, anDof, anOff;
31270034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
313cb1e1211SMatthew G Knepley 
314fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
315cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
316cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
317a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
31876185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
31976185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
320cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
321cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
322cb1e1211SMatthew G Knepley 
323cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(leafSectionAdj, d, &aoff);CHKERRQ(ierr);
324cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
325cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
326cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
327cb1e1211SMatthew G Knepley 
328cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
329cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
330cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
331cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
332cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
333cb1e1211SMatthew G Knepley           adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd;
334cb1e1211SMatthew G Knepley           ++i;
335cb1e1211SMatthew G Knepley         }
336cb1e1211SMatthew G Knepley       }
33776185916SToby Isaac       for (q = 0; q < anDof; q++) {
33876185916SToby Isaac         adj[aoff+i] = anchorAdj[anOff+q];
33976185916SToby Isaac         ++i;
34076185916SToby Isaac       }
341cb1e1211SMatthew G Knepley     }
342cb1e1211SMatthew G Knepley   }
343cb1e1211SMatthew G Knepley   /* Debugging */
344cb1e1211SMatthew G Knepley   if (debug) {
345cb1e1211SMatthew G Knepley     IS tmp;
346cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Leaf adjacency indices\n");CHKERRQ(ierr);
347cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
348cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
3490a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
350cb1e1211SMatthew G Knepley   }
351543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
352cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(rootSectionAdj, &adjSize);CHKERRQ(ierr);
353785e854fSJed Brown   ierr = PetscMalloc1(adjSize, &rootAdj);CHKERRQ(ierr);
354cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
35547a6104aSMatthew G. Knepley   if (doComm) {
356543482b8SMatthew G. Knepley     const PetscInt *indegree;
357543482b8SMatthew G. Knepley     PetscInt       *remoteadj, radjsize = 0;
358543482b8SMatthew G. Knepley 
359543482b8SMatthew G. Knepley     ierr = PetscSFComputeDegreeBegin(sfAdj, &indegree);CHKERRQ(ierr);
360543482b8SMatthew G. Knepley     ierr = PetscSFComputeDegreeEnd(sfAdj, &indegree);CHKERRQ(ierr);
361543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
362543482b8SMatthew G. Knepley     ierr = PetscMalloc1(radjsize, &remoteadj);CHKERRQ(ierr);
363543482b8SMatthew G. Knepley     ierr = PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr);
364543482b8SMatthew G. Knepley     ierr = PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr);
3656dba2905SMatthew G. Knepley     for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p-1])) {
366543482b8SMatthew G. Knepley       PetscInt s;
3676dba2905SMatthew G. Knepley       for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l+s] = remoteadj[r];
368543482b8SMatthew G. Knepley     }
369543482b8SMatthew G. Knepley     if (r != radjsize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", r, radjsize);
3706dba2905SMatthew G. Knepley     if (l != adjSize)  SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", l, adjSize);
371543482b8SMatthew G. Knepley     ierr = PetscFree(remoteadj);CHKERRQ(ierr);
372cb1e1211SMatthew G Knepley   }
373cb1e1211SMatthew G Knepley   ierr = PetscSFDestroy(&sfAdj);CHKERRQ(ierr);
374cb1e1211SMatthew G Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
375cb1e1211SMatthew G Knepley   /* Debugging */
376cb1e1211SMatthew G Knepley   if (debug) {
377cb1e1211SMatthew G Knepley     IS tmp;
378cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after gather\n");CHKERRQ(ierr);
379cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
380cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
3810a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
382cb1e1211SMatthew G Knepley   }
383cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
384cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
38576185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
386cb1e1211SMatthew G Knepley 
387cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
388cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
389cb1e1211SMatthew G Knepley     if (!dof) continue;
390cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
391cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
392a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
39376185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
39476185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
395cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
396cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
397cb1e1211SMatthew G Knepley 
398cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
399cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
400cb1e1211SMatthew G Knepley       i    = adof-1;
40176185916SToby Isaac       for (q = 0; q < anDof; q++) {
40276185916SToby Isaac         rootAdj[aoff+i] = anchorAdj[anOff+q];
40376185916SToby Isaac         --i;
40476185916SToby Isaac       }
405cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
406cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
407cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
408cb1e1211SMatthew G Knepley 
409cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
410cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
411cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
412cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
413cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
414cb1e1211SMatthew G Knepley           rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
415cb1e1211SMatthew G Knepley           --i;
416cb1e1211SMatthew G Knepley         }
417cb1e1211SMatthew G Knepley       }
418cb1e1211SMatthew G Knepley     }
419cb1e1211SMatthew G Knepley   }
420cb1e1211SMatthew G Knepley   /* Debugging */
421cb1e1211SMatthew G Knepley   if (debug) {
422cb1e1211SMatthew G Knepley     IS tmp;
423cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices\n");CHKERRQ(ierr);
424cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
425cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
4260a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
427cb1e1211SMatthew G Knepley   }
428cb1e1211SMatthew G Knepley   /* Compress indices */
429cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
430cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
431cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
432cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
433cb1e1211SMatthew G Knepley 
434cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
435cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
436cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
437cb1e1211SMatthew G Knepley     if (!dof) continue;
438cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
439cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
440cb1e1211SMatthew G Knepley     for (d = off; d < off+dof-cdof; ++d) {
441cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
442cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
443cb1e1211SMatthew G Knepley       ierr = PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]);CHKERRQ(ierr);
444cb1e1211SMatthew G Knepley       ierr = PetscSectionSetDof(rootSectionAdj, d, adof);CHKERRQ(ierr);
445cb1e1211SMatthew G Knepley     }
446cb1e1211SMatthew G Knepley   }
447cb1e1211SMatthew G Knepley   /* Debugging */
448cb1e1211SMatthew G Knepley   if (debug) {
449cb1e1211SMatthew G Knepley     IS tmp;
450cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n");CHKERRQ(ierr);
451cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
452cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after compression\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   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
458cb1e1211SMatthew G Knepley   ierr = PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd);CHKERRQ(ierr);
459cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &sectionAdj);CHKERRQ(ierr);
460cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd);CHKERRQ(ierr);
461cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
46276185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
463cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
464cb1e1211SMatthew G Knepley 
465cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
466cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
467cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
468cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
469cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
470cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
471cb1e1211SMatthew G Knepley 
472cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
473cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
474cb1e1211SMatthew G Knepley       if (ldof > 0) {
475cb1e1211SMatthew G Knepley         /* We do not own this point */
476cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
477cb1e1211SMatthew G Knepley         ierr = PetscSectionSetDof(sectionAdj, goff+d, rdof);CHKERRQ(ierr);
478cb1e1211SMatthew G Knepley       } else {
479cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
480cb1e1211SMatthew G Knepley       }
481cb1e1211SMatthew G Knepley     }
482cb1e1211SMatthew G Knepley     if (found) continue;
483cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
484cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
485a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
486cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
487cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
488cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof, noff;
489cb1e1211SMatthew G Knepley 
490cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
491cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
492cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
493cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(section, padj, &noff);CHKERRQ(ierr);
494cb1e1211SMatthew G Knepley       for (d = goff; d < goff+dof-cdof; ++d) {
495cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(sectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
496cb1e1211SMatthew G Knepley       }
497cb1e1211SMatthew G Knepley     }
49876185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
49976185916SToby Isaac     if (anDof) {
50076185916SToby Isaac       for (d = goff; d < goff+dof-cdof; ++d) {
50176185916SToby Isaac         ierr = PetscSectionAddDof(sectionAdj, d, anDof);CHKERRQ(ierr);
50276185916SToby Isaac       }
50376185916SToby Isaac     }
504cb1e1211SMatthew G Knepley   }
505cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(sectionAdj);CHKERRQ(ierr);
506cb1e1211SMatthew G Knepley   if (debug) {
507cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation:\n");CHKERRQ(ierr);
508cb1e1211SMatthew G Knepley     ierr = PetscSectionView(sectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
509cb1e1211SMatthew G Knepley   }
510cb1e1211SMatthew G Knepley   /* Get adjacent indices */
511cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(sectionAdj, &numCols);CHKERRQ(ierr);
512785e854fSJed Brown   ierr = PetscMalloc1(numCols, &cols);CHKERRQ(ierr);
513cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
51476185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
515cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
516cb1e1211SMatthew G Knepley 
517cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
518cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
519cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
520cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
521cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
522cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
523cb1e1211SMatthew G Knepley 
524cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
525cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
526cb1e1211SMatthew G Knepley       if (ldof > 0) {
527cb1e1211SMatthew G Knepley         /* We do not own this point */
528cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
529cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
530cb1e1211SMatthew G Knepley 
531cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionAdj, goff+d, &aoff);CHKERRQ(ierr);
532cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(rootSectionAdj, off+d, &roff);CHKERRQ(ierr);
533cb1e1211SMatthew G Knepley         ierr = PetscMemcpy(&cols[aoff], &rootAdj[roff], rdof * sizeof(PetscInt));CHKERRQ(ierr);
534cb1e1211SMatthew G Knepley       } else {
535cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
536cb1e1211SMatthew G Knepley       }
537cb1e1211SMatthew G Knepley     }
538cb1e1211SMatthew G Knepley     if (found) continue;
539a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
54076185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
54176185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
542cb1e1211SMatthew G Knepley     for (d = goff; d < goff+dof-cdof; ++d) {
543cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
544cb1e1211SMatthew G Knepley 
545cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(sectionAdj, d, &adof);CHKERRQ(ierr);
546cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(sectionAdj, d, &aoff);CHKERRQ(ierr);
547cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
548cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
549cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
550cb1e1211SMatthew G Knepley         const PetscInt *ncind;
551cb1e1211SMatthew G Knepley 
552cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
553cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
554cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
555cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
556cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintIndices(section, padj, &ncind);CHKERRQ(ierr);
557cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
558cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd, ++i) {
559cb1e1211SMatthew G Knepley           cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
560cb1e1211SMatthew G Knepley         }
561cb1e1211SMatthew G Knepley       }
56276185916SToby Isaac       for (q = 0; q < anDof; q++, i++) {
56376185916SToby Isaac         cols[aoff+i] = anchorAdj[anOff + q];
56476185916SToby Isaac       }
565cb1e1211SMatthew 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);
566cb1e1211SMatthew G Knepley     }
567cb1e1211SMatthew G Knepley   }
56876185916SToby Isaac   ierr = PetscSectionDestroy(&anchorSectionAdj);CHKERRQ(ierr);
569cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&leafSectionAdj);CHKERRQ(ierr);
570cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&rootSectionAdj);CHKERRQ(ierr);
57176185916SToby Isaac   ierr = PetscFree(anchorAdj);CHKERRQ(ierr);
572cb1e1211SMatthew G Knepley   ierr = PetscFree(rootAdj);CHKERRQ(ierr);
57370034214SMatthew G. Knepley   ierr = PetscFree(tmpAdj);CHKERRQ(ierr);
574cb1e1211SMatthew G Knepley   /* Debugging */
575cb1e1211SMatthew G Knepley   if (debug) {
576cb1e1211SMatthew G Knepley     IS tmp;
577cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Column indices\n");CHKERRQ(ierr);
578cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
579cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
5800a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
581cb1e1211SMatthew G Knepley   }
582a9fb6443SMatthew G. Knepley 
583a9fb6443SMatthew G. Knepley   *sA     = sectionAdj;
584a9fb6443SMatthew G. Knepley   *colIdx = cols;
585a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
586a9fb6443SMatthew G. Knepley }
587a9fb6443SMatthew G. Knepley 
588a9fb6443SMatthew G. Knepley #undef __FUNCT__
589a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexUpdateAllocation_Static"
590e101f074SMatthew G. Knepley PetscErrorCode DMPlexUpdateAllocation_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[])
591a9fb6443SMatthew G. Knepley {
592e101f074SMatthew G. Knepley   PetscSection   section;
593a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p;
594a9fb6443SMatthew G. Knepley   PetscErrorCode ierr;
595a9fb6443SMatthew G. Knepley 
596a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
597a9fb6443SMatthew G. Knepley   /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */
598cb1e1211SMatthew G Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
599a9fb6443SMatthew G. Knepley   if (rStart%bs || rEnd%bs) SETERRQ3(PetscObjectComm((PetscObject) rLayout), PETSC_ERR_ARG_WRONG, "Invalid layout [%d, %d) for matrix, must be divisible by block size %d", rStart, rEnd, bs);
600a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
601e101f074SMatthew G. Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
602a9fb6443SMatthew G. Knepley     ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
603a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
604294b7009SMatthew G. Knepley       PetscInt rS, rE;
605a9fb6443SMatthew G. Knepley 
606294b7009SMatthew G. Knepley       ierr = DMPlexGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE);CHKERRQ(ierr);
607a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
608a9fb6443SMatthew G. Knepley         PetscInt numCols, cStart, c;
609a9fb6443SMatthew G. Knepley 
610a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
611a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
612a9fb6443SMatthew G. Knepley         for (c = cStart; c < cStart+numCols; ++c) {
613a9fb6443SMatthew G. Knepley           if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
614a9fb6443SMatthew G. Knepley             ++dnz[r-rStart];
615a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++dnzu[r-rStart];
616a9fb6443SMatthew G. Knepley           } else {
617a9fb6443SMatthew G. Knepley             ++onz[r-rStart];
618a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++onzu[r-rStart];
619a9fb6443SMatthew G. Knepley           }
620a9fb6443SMatthew G. Knepley         }
621a9fb6443SMatthew G. Knepley       }
622a9fb6443SMatthew G. Knepley     }
623a9fb6443SMatthew G. Knepley   } else {
624cb1e1211SMatthew G Knepley     /* Only loop over blocks of rows */
625cb1e1211SMatthew G Knepley     for (r = rStart/bs; r < rEnd/bs; ++r) {
626cb1e1211SMatthew G Knepley       const PetscInt row = r*bs;
627cb1e1211SMatthew G Knepley       PetscInt       numCols, cStart, c;
628cb1e1211SMatthew G Knepley 
629cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(sectionAdj, row, &numCols);CHKERRQ(ierr);
630cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(sectionAdj, row, &cStart);CHKERRQ(ierr);
631cb1e1211SMatthew G Knepley       for (c = cStart; c < cStart+numCols; ++c) {
632*e7bcfa36SMatthew G. Knepley         if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
633*e7bcfa36SMatthew G. Knepley           ++dnz[r-rStart/bs];
634*e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++dnzu[r-rStart/bs];
635cb1e1211SMatthew G Knepley         } else {
636*e7bcfa36SMatthew G. Knepley           ++onz[r-rStart/bs];
637*e7bcfa36SMatthew G. Knepley           if (cols[c] >= row) ++onzu[r-rStart/bs];
638cb1e1211SMatthew G Knepley         }
639cb1e1211SMatthew G Knepley       }
640cb1e1211SMatthew G Knepley     }
641a9fb6443SMatthew G. Knepley     for (r = 0; r < (rEnd - rStart)/bs; ++r) {
642cb1e1211SMatthew G Knepley       dnz[r]  /= bs;
643cb1e1211SMatthew G Knepley       onz[r]  /= bs;
644cb1e1211SMatthew G Knepley       dnzu[r] /= bs;
645cb1e1211SMatthew G Knepley       onzu[r] /= bs;
646cb1e1211SMatthew G Knepley     }
647cb1e1211SMatthew G Knepley   }
648a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
649a9fb6443SMatthew G. Knepley }
650a9fb6443SMatthew G. Knepley 
651a9fb6443SMatthew G. Knepley #undef __FUNCT__
652a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexFillMatrix_Static"
653e101f074SMatthew G. Knepley PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A)
654a9fb6443SMatthew G. Knepley {
655e101f074SMatthew G. Knepley   PetscSection   section;
656a9fb6443SMatthew G. Knepley   PetscScalar   *values;
657a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0;
658a9fb6443SMatthew G. Knepley   PetscErrorCode ierr;
659a9fb6443SMatthew G. Knepley 
660a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
661a9fb6443SMatthew G. Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
662a9fb6443SMatthew G. Knepley   for (r = rStart; r < rEnd; ++r) {
663a9fb6443SMatthew G. Knepley     ierr      = PetscSectionGetDof(sectionAdj, r, &len);CHKERRQ(ierr);
664a9fb6443SMatthew G. Knepley     maxRowLen = PetscMax(maxRowLen, len);
665a9fb6443SMatthew G. Knepley   }
666a9fb6443SMatthew G. Knepley   ierr = PetscCalloc1(maxRowLen, &values);CHKERRQ(ierr);
667a9fb6443SMatthew G. Knepley   if (f >=0 && bs == 1) {
668e101f074SMatthew G. Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
669a9fb6443SMatthew G. Knepley     ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
670a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
671294b7009SMatthew G. Knepley       PetscInt rS, rE;
672a9fb6443SMatthew G. Knepley 
673294b7009SMatthew G. Knepley       ierr = DMPlexGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE);CHKERRQ(ierr);
674a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
6757e01ee02SMatthew G. Knepley         PetscInt numCols, cStart;
676a9fb6443SMatthew G. Knepley 
677a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
678a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
679a9fb6443SMatthew G. Knepley         ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr);
680a9fb6443SMatthew G. Knepley       }
681a9fb6443SMatthew G. Knepley     }
682a9fb6443SMatthew G. Knepley   } else {
683a9fb6443SMatthew G. Knepley     for (r = rStart; r < rEnd; ++r) {
684a9fb6443SMatthew G. Knepley       PetscInt numCols, cStart;
685a9fb6443SMatthew G. Knepley 
686a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
687a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
688a9fb6443SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr);
689a9fb6443SMatthew G. Knepley     }
690a9fb6443SMatthew G. Knepley   }
691a9fb6443SMatthew G. Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
692a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
693a9fb6443SMatthew G. Knepley }
694a9fb6443SMatthew G. Knepley 
695a9fb6443SMatthew G. Knepley #undef __FUNCT__
696a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexPreallocateOperator"
697e101f074SMatthew G. Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
698a9fb6443SMatthew G. Knepley {
699a9fb6443SMatthew G. Knepley   MPI_Comm       comm;
700a9fb6443SMatthew G. Knepley   PetscDS        prob;
701a9fb6443SMatthew G. Knepley   MatType        mtype;
702a9fb6443SMatthew G. Knepley   PetscSF        sf, sfDof;
703e101f074SMatthew G. Knepley   PetscSection   section;
704a9fb6443SMatthew G. Knepley   PetscInt      *remoteOffsets;
705a9fb6443SMatthew G. Knepley   PetscSection   sectionAdj[4] = {NULL, NULL, NULL, NULL};
706a9fb6443SMatthew G. Knepley   PetscInt      *cols[4]       = {NULL, NULL, NULL, NULL};
707a9fb6443SMatthew G. Knepley   PetscBool      useCone, useClosure;
7087e01ee02SMatthew G. Knepley   PetscInt       Nf, f, idx, locRows;
709a9fb6443SMatthew G. Knepley   PetscLayout    rLayout;
710a9fb6443SMatthew G. Knepley   PetscBool      isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE;
711a9fb6443SMatthew G. Knepley   PetscErrorCode ierr;
712a9fb6443SMatthew G. Knepley 
713a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
714a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
715a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(A, MAT_CLASSID, 9);
716a9fb6443SMatthew G. Knepley   if (dnz)  PetscValidPointer(dnz,5);  if (onz)  PetscValidPointer(onz,6);
717a9fb6443SMatthew G. Knepley   if (dnzu) PetscValidPointer(dnzu,7); if (onzu) PetscValidPointer(onzu,8);
718a9fb6443SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
719a9fb6443SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
720e101f074SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
721c5929fdfSBarry Smith   ierr = PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr);
722a9fb6443SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
723a9fb6443SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr);
724a9fb6443SMatthew G. Knepley   /* Create dof SF based on point SF */
725a9fb6443SMatthew G. Knepley   if (debug) {
726e101f074SMatthew G. Knepley     PetscSection section, sectionGlobal;
727a9fb6443SMatthew G. Knepley     PetscSF      sf;
728a9fb6443SMatthew G. Knepley 
729a9fb6443SMatthew G. Knepley     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
730e101f074SMatthew G. Knepley     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
731e101f074SMatthew G. Knepley     ierr = DMGetDefaultGlobalSection(dm, &sectionGlobal);CHKERRQ(ierr);
732a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Input Section for Preallocation:\n");CHKERRQ(ierr);
733a9fb6443SMatthew G. Knepley     ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
734a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Input Global Section for Preallocation:\n");CHKERRQ(ierr);
735a9fb6443SMatthew G. Knepley     ierr = PetscSectionView(sectionGlobal, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
736a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Input SF for Preallocation:\n");CHKERRQ(ierr);
737a9fb6443SMatthew G. Knepley     ierr = PetscSFView(sf, NULL);CHKERRQ(ierr);
738a9fb6443SMatthew G. Knepley   }
739a9fb6443SMatthew G. Knepley   ierr = PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets);CHKERRQ(ierr);
740a9fb6443SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof);CHKERRQ(ierr);
7418a713f3bSLawrence Mitchell   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
742a9fb6443SMatthew G. Knepley   if (debug) {
743a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Dof SF for Preallocation:\n");CHKERRQ(ierr);
744a9fb6443SMatthew G. Knepley     ierr = PetscSFView(sfDof, NULL);CHKERRQ(ierr);
745a9fb6443SMatthew G. Knepley   }
746a9fb6443SMatthew G. Knepley   /* Create allocation vectors from adjacency graph */
747a9fb6443SMatthew G. Knepley   ierr = MatGetLocalSize(A, &locRows, NULL);CHKERRQ(ierr);
748a9fb6443SMatthew G. Knepley   ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr);
749a9fb6443SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
750a9fb6443SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
751a9fb6443SMatthew G. Knepley   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
752a9fb6443SMatthew G. Knepley   /* There are 4 types of adjacency */
753a9fb6443SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
754a9fb6443SMatthew G. Knepley   if (Nf < 1 || bs > 1) {
755a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
756a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
757a9fb6443SMatthew G. Knepley     idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
758e101f074SMatthew G. Knepley     ierr = DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]);CHKERRQ(ierr);
759e101f074SMatthew G. Knepley     ierr = DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu);CHKERRQ(ierr);
760a9fb6443SMatthew G. Knepley   } else {
761a9fb6443SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
762a9fb6443SMatthew G. Knepley       ierr = PetscDSGetAdjacency(prob, f, &useCone, &useClosure);CHKERRQ(ierr);
763a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
764e101f074SMatthew G. Knepley       if (!sectionAdj[idx]) {ierr = DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]);CHKERRQ(ierr);}
765e101f074SMatthew G. Knepley       ierr = DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu);CHKERRQ(ierr);
766a9fb6443SMatthew G. Knepley     }
767a9fb6443SMatthew G. Knepley   }
768a9fb6443SMatthew G. Knepley   ierr = PetscSFDestroy(&sfDof);CHKERRQ(ierr);
769cb1e1211SMatthew G Knepley   /* Set matrix pattern */
770cb1e1211SMatthew G Knepley   ierr = MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu);CHKERRQ(ierr);
771cb1e1211SMatthew G Knepley   ierr = MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
77289545effSMatthew G. Knepley   /* Check for symmetric storage */
77389545effSMatthew G. Knepley   ierr = MatGetType(A, &mtype);CHKERRQ(ierr);
77489545effSMatthew G. Knepley   ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr);
77589545effSMatthew G. Knepley   ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr);
77689545effSMatthew G. Knepley   ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr);
77789545effSMatthew G. Knepley   if (isSymBlock || isSymSeqBlock || isSymMPIBlock) {ierr = MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);}
778cb1e1211SMatthew G Knepley   /* Fill matrix with zeros */
779cb1e1211SMatthew G Knepley   if (fillMatrix) {
780a9fb6443SMatthew G. Knepley     if (Nf < 1 || bs > 1) {
781a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
782a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
783a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
784e101f074SMatthew G. Knepley       ierr = DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A);CHKERRQ(ierr);
785a9fb6443SMatthew G. Knepley     } else {
786a9fb6443SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
787a9fb6443SMatthew G. Knepley         ierr = PetscDSGetAdjacency(prob, f, &useCone, &useClosure);CHKERRQ(ierr);
788a9fb6443SMatthew G. Knepley         idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
789e101f074SMatthew G. Knepley         ierr = DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A);CHKERRQ(ierr);
790cb1e1211SMatthew G Knepley       }
791cb1e1211SMatthew G Knepley     }
792cb1e1211SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
793cb1e1211SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
794cb1e1211SMatthew G Knepley   }
795a9fb6443SMatthew G. Knepley   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
796a9fb6443SMatthew G. Knepley   for (idx = 0; idx < 4; ++idx) {ierr = PetscSectionDestroy(&sectionAdj[idx]);CHKERRQ(ierr); ierr = PetscFree(cols[idx]);CHKERRQ(ierr);}
797a72f3261SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr);
798cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
799cb1e1211SMatthew G Knepley }
800cb1e1211SMatthew G Knepley 
801cb1e1211SMatthew G Knepley #if 0
802cb1e1211SMatthew G Knepley #undef __FUNCT__
803cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator_2"
804cb1e1211SMatthew 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)
805cb1e1211SMatthew G Knepley {
806cb1e1211SMatthew G Knepley   PetscInt       *tmpClosure,*tmpAdj,*visits;
807cb1e1211SMatthew G Knepley   PetscInt        c,cStart,cEnd,pStart,pEnd;
808cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
809cb1e1211SMatthew G Knepley 
810cb1e1211SMatthew G Knepley   PetscFunctionBegin;
811c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
812cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
813cb1e1211SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
814cb1e1211SMatthew G Knepley 
815e7b80042SMatthew G. Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1));
816cb1e1211SMatthew G Knepley 
817cb1e1211SMatthew G Knepley   ierr    = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
818cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
819cb1e1211SMatthew G Knepley 
820dcca6d9dSJed Brown   ierr = PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits);CHKERRQ(ierr);
821cb1e1211SMatthew G Knepley   ierr = PetscMemzero(lvisits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr);
822cb1e1211SMatthew G Knepley   ierr = PetscMemzero(visits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr);
823cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
824cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
825cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
826cb1e1211SMatthew G Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support);CHKERRQ(ierr);
827cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
828cb1e1211SMatthew G Knepley   }
829cb1e1211SMatthew G Knepley   ierr = PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr);
830cb1e1211SMatthew G Knepley   ierr = PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr);
831cb1e1211SMatthew G Knepley   ierr = PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr);
832cb1e1211SMatthew G Knepley   ierr = PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr);
833cb1e1211SMatthew G Knepley 
834cb1e1211SMatthew G Knepley   ierr = PetscSFGetRanks();CHKERRQ(ierr);
835cb1e1211SMatthew G Knepley 
836cb1e1211SMatthew G Knepley 
837dcca6d9dSJed Brown   ierr = PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner);CHKERRQ(ierr);
838cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
839cb1e1211SMatthew G Knepley     ierr = PetscMemzero(cellmat,maxClosureSize*maxClosureSize*sizeof(PetscInt));CHKERRQ(ierr);
840cb1e1211SMatthew G Knepley     /*
841cb1e1211SMatthew G Knepley      Depth-first walk of transitive closure.
842cb1e1211SMatthew 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.
843cb1e1211SMatthew G Knepley      This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
844cb1e1211SMatthew G Knepley      */
845cb1e1211SMatthew G Knepley   }
846cb1e1211SMatthew G Knepley 
847cb1e1211SMatthew G Knepley   ierr = PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM);CHKERRQ(ierr);
848cb1e1211SMatthew G Knepley   ierr = PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM);CHKERRQ(ierr);
849cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
850cb1e1211SMatthew G Knepley }
851cb1e1211SMatthew G Knepley #endif
852