xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision a9fb644384c31e1fad7fabb80d80dc0b100418b9)
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>
4*a9fb6443SMatthew 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() */
9*a9fb6443SMatthew G. Knepley static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscSection section, PetscSection sectionGlobal, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[])
1076185916SToby Isaac {
1176185916SToby Isaac   PetscInt       pStart, pEnd;
1276185916SToby Isaac   PetscSection   adjSec, aSec;
1376185916SToby Isaac   IS             aIS;
1476185916SToby Isaac   PetscErrorCode ierr;
1576185916SToby Isaac 
1676185916SToby Isaac   PetscFunctionBegin;
1776185916SToby Isaac 
1876185916SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&adjSec);CHKERRQ(ierr);
1976185916SToby Isaac   ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr);
2076185916SToby Isaac   ierr = PetscSectionSetChart(adjSec,pStart,pEnd);CHKERRQ(ierr);
2176185916SToby Isaac 
22a17985deSToby Isaac   ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
2376185916SToby Isaac   if (aSec) {
2476185916SToby Isaac     const PetscInt *anchors;
2576185916SToby Isaac     PetscInt       p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize;
2676185916SToby Isaac     PetscInt       *tmpAdjP = NULL, *tmpAdjQ = NULL;
2776185916SToby Isaac     PetscSection   inverseSec;
2876185916SToby Isaac 
2976185916SToby Isaac     /* invert the constraint-to-anchor map */
3076185916SToby Isaac     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)aSec),&inverseSec);CHKERRQ(ierr);
3176185916SToby Isaac     ierr = PetscSectionSetChart(inverseSec,pStart,pEnd);CHKERRQ(ierr);
3276185916SToby Isaac     ierr = ISGetLocalSize(aIS, &aSize);CHKERRQ(ierr);
3376185916SToby Isaac     ierr = ISGetIndices(aIS, &anchors);CHKERRQ(ierr);
3476185916SToby Isaac 
3576185916SToby Isaac     for (p = 0; p < aSize; p++) {
3676185916SToby Isaac       PetscInt a = anchors[p];
3776185916SToby Isaac 
3876185916SToby Isaac       ierr = PetscSectionAddDof(inverseSec,a,1);CHKERRQ(ierr);
3976185916SToby Isaac     }
4076185916SToby Isaac     ierr = PetscSectionSetUp(inverseSec);CHKERRQ(ierr);
4176185916SToby Isaac     ierr = PetscSectionGetStorageSize(inverseSec,&iSize);CHKERRQ(ierr);
4276185916SToby Isaac     ierr = PetscMalloc1(iSize,&inverse);CHKERRQ(ierr);
4376185916SToby Isaac     ierr = PetscCalloc1(pEnd-pStart,&offsets);CHKERRQ(ierr);
4476185916SToby Isaac     ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
4576185916SToby Isaac     for (p = aStart; p < aEnd; p++) {
4676185916SToby Isaac       PetscInt dof, off;
4776185916SToby Isaac 
4876185916SToby Isaac       ierr = PetscSectionGetDof(aSec, p, &dof);CHKERRQ(ierr);
4976185916SToby Isaac       ierr = PetscSectionGetOffset(aSec, p, &off);CHKERRQ(ierr);
5076185916SToby Isaac 
5176185916SToby Isaac       for (q = 0; q < dof; q++) {
5276185916SToby Isaac         PetscInt iOff;
5376185916SToby Isaac 
5476185916SToby Isaac         a = anchors[off + q];
5576185916SToby Isaac         ierr = PetscSectionGetOffset(inverseSec, a, &iOff);CHKERRQ(ierr);
5676185916SToby Isaac         inverse[iOff + offsets[a-pStart]++] = p;
5776185916SToby Isaac       }
5876185916SToby Isaac     }
5976185916SToby Isaac     ierr = ISRestoreIndices(aIS, &anchors);CHKERRQ(ierr);
6076185916SToby Isaac     ierr = PetscFree(offsets);CHKERRQ(ierr);
6176185916SToby Isaac 
6276185916SToby Isaac     /* construct anchorAdj and adjSec
6376185916SToby Isaac      *
6476185916SToby Isaac      * loop over anchors:
6576185916SToby Isaac      *   construct anchor adjacency
6676185916SToby Isaac      *   loop over constrained:
6776185916SToby Isaac      *     construct constrained adjacency
6876185916SToby Isaac      *     if not in anchor adjacency, add to dofs
6976185916SToby Isaac      * setup adjSec, allocate anchorAdj
7076185916SToby Isaac      * loop over anchors:
7176185916SToby Isaac      *   construct anchor adjacency
7276185916SToby Isaac      *   loop over constrained:
7376185916SToby Isaac      *     construct constrained adjacency
7476185916SToby Isaac      *     if not in anchor adjacency
7576185916SToby Isaac      *       if not already in list, put in list
7676185916SToby Isaac      *   sort, unique, reduce dof count
7776185916SToby Isaac      * optional: compactify
7876185916SToby Isaac      */
7976185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
8076185916SToby Isaac       PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE;
8176185916SToby Isaac 
8276185916SToby Isaac       ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr);
8376185916SToby Isaac       if (!iDof) continue;
8476185916SToby Isaac       ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr);
85*a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP);CHKERRQ(ierr);
8676185916SToby Isaac       for (i = 0; i < iDof; i++) {
8776185916SToby Isaac         PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE;
8876185916SToby Isaac 
8976185916SToby Isaac         q = inverse[iOff + i];
90*a9fb6443SMatthew G. Knepley         ierr = DMPlexGetAdjacency_Internal(dm,q,useCone,useClosure,PETSC_TRUE,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr);
9176185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
9276185916SToby Isaac           qAdj = tmpAdjQ[r];
9376185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
9476185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
9576185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
9676185916SToby Isaac           }
9776185916SToby Isaac           if (s < numAdjP) continue;
9876185916SToby Isaac           ierr  = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr);
9976185916SToby Isaac           ierr  = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr);
10076185916SToby Isaac           iNew += qAdjDof - qAdjCDof;
10176185916SToby Isaac         }
10276185916SToby Isaac         ierr = PetscSectionAddDof(adjSec,p,iNew);CHKERRQ(ierr);
10376185916SToby Isaac       }
10476185916SToby Isaac     }
10576185916SToby Isaac 
10676185916SToby Isaac     ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr);
10776185916SToby Isaac     ierr = PetscSectionGetStorageSize(adjSec,&adjSize);CHKERRQ(ierr);
10876185916SToby Isaac     ierr = PetscMalloc1(adjSize,&adj);CHKERRQ(ierr);
10976185916SToby Isaac 
11076185916SToby Isaac     for (p = pStart; p < pEnd; p++) {
11176185916SToby Isaac       PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE;
11276185916SToby Isaac 
11376185916SToby Isaac       ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr);
11476185916SToby Isaac       if (!iDof) continue;
11576185916SToby Isaac       ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr);
116*a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP);CHKERRQ(ierr);
11776185916SToby Isaac       ierr = PetscSectionGetDof(adjSec,p,&aDof);CHKERRQ(ierr);
11876185916SToby Isaac       ierr = PetscSectionGetOffset(adjSec,p,&aOff);CHKERRQ(ierr);
11976185916SToby Isaac       aOffOrig = aOff;
12076185916SToby Isaac       for (i = 0; i < iDof; i++) {
12176185916SToby Isaac         PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE;
12276185916SToby Isaac 
12376185916SToby Isaac         q = inverse[iOff + i];
124*a9fb6443SMatthew G. Knepley         ierr = DMPlexGetAdjacency_Internal(dm,q,useCone,useClosure,PETSC_TRUE,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr);
12576185916SToby Isaac         for (r = 0; r < numAdjQ; r++) {
12676185916SToby Isaac           qAdj = tmpAdjQ[r];
12776185916SToby Isaac           if ((qAdj < pStart) || (qAdj >= pEnd)) continue;
12876185916SToby Isaac           for (s = 0; s < numAdjP; s++) {
12976185916SToby Isaac             if (qAdj == tmpAdjP[s]) break;
13076185916SToby Isaac           }
13176185916SToby Isaac           if (s < numAdjP) continue;
13276185916SToby Isaac           ierr  = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr);
13376185916SToby Isaac           ierr  = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr);
13476185916SToby Isaac           ierr  = PetscSectionGetOffset(sectionGlobal,qAdj,&qAdjOff);CHKERRQ(ierr);
13576185916SToby Isaac           for (nd = 0; nd < qAdjDof-qAdjCDof; ++nd) {
13676185916SToby Isaac             adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff+1) : qAdjOff) + nd;
13776185916SToby Isaac           }
13876185916SToby Isaac         }
13976185916SToby Isaac       }
14076185916SToby Isaac       ierr = PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]);CHKERRQ(ierr);
14176185916SToby Isaac       ierr = PetscSectionSetDof(adjSec,p,aDof);CHKERRQ(ierr);
14276185916SToby Isaac     }
14376185916SToby Isaac     *anchorAdj = adj;
14476185916SToby Isaac 
14576185916SToby Isaac     /* clean up */
14676185916SToby Isaac     ierr = PetscSectionDestroy(&inverseSec);CHKERRQ(ierr);
14776185916SToby Isaac     ierr = PetscFree(inverse);CHKERRQ(ierr);
14876185916SToby Isaac     ierr = PetscFree(tmpAdjP);CHKERRQ(ierr);
14976185916SToby Isaac     ierr = PetscFree(tmpAdjQ);CHKERRQ(ierr);
15076185916SToby Isaac   }
15176185916SToby Isaac   else {
15276185916SToby Isaac     *anchorAdj = NULL;
15376185916SToby Isaac     ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr);
15476185916SToby Isaac   }
15576185916SToby Isaac   *anchorSectionAdj = adjSec;
15676185916SToby Isaac   PetscFunctionReturn(0);
15776185916SToby Isaac }
15876185916SToby Isaac 
15976185916SToby Isaac #undef __FUNCT__
160*a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateAdjacencySection_Static"
161*a9fb6443SMatthew G. Knepley PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscSection section, PetscSection sectionGlobal, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx)
162cb1e1211SMatthew G Knepley {
163cb1e1211SMatthew G Knepley   MPI_Comm           comm;
164*a9fb6443SMatthew G. Knepley   PetscMPIInt        size;
165*a9fb6443SMatthew G. Knepley   PetscBool          doCommLocal, doComm, debug = PETSC_FALSE;
166*a9fb6443SMatthew G. Knepley   PetscSF            sf, sfAdj;
16776185916SToby Isaac   PetscSection       leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj;
168*a9fb6443SMatthew G. Knepley   PetscInt           nroots, nleaves, l, p, r;
169cb1e1211SMatthew G Knepley   const PetscInt    *leaves;
170cb1e1211SMatthew G Knepley   const PetscSFNode *remotes;
171cb1e1211SMatthew G Knepley   PetscInt           dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols;
1728821704fSMatthew G. Knepley   PetscInt          *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets;
17370034214SMatthew G. Knepley   PetscInt           adjSize;
174cb1e1211SMatthew G Knepley   PetscErrorCode     ierr;
175cb1e1211SMatthew G Knepley 
176cb1e1211SMatthew G Knepley   PetscFunctionBegin;
177cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
178cb1e1211SMatthew G Knepley   ierr = PetscOptionsGetBool(NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr);
179cb1e1211SMatthew G Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
180c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
181cb1e1211SMatthew G Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
18247a6104aSMatthew G. Knepley   ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
18347a6104aSMatthew G. Knepley   doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE;
18447a6104aSMatthew G. Knepley   ierr = MPI_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr);
185cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
186cb1e1211SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
187cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &numDof);CHKERRQ(ierr);
188cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &leafSectionAdj);CHKERRQ(ierr);
189cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(leafSectionAdj, 0, numDof);CHKERRQ(ierr);
190cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &rootSectionAdj);CHKERRQ(ierr);
191cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(rootSectionAdj, 0, numDof);CHKERRQ(ierr);
192cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
193cb1e1211SMatthew G Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes);CHKERRQ(ierr);
194cb1e1211SMatthew G Knepley   /*
195964bf7afSMatthew G. Knepley    section        - maps points to (# dofs, local dofs)
196964bf7afSMatthew G. Knepley    sectionGlobal  - maps points to (# dofs, global dofs)
197964bf7afSMatthew G. Knepley    leafSectionAdj - maps unowned local dofs to # adj dofs
198964bf7afSMatthew G. Knepley    rootSectionAdj - maps   owned local dofs to # adj dofs
199964bf7afSMatthew G. Knepley    adj            - adj global dofs indexed by leafSectionAdj
200964bf7afSMatthew G. Knepley    rootAdj        - adj global dofs indexed by rootSectionAdj
201964bf7afSMatthew G. Knepley    sf    - describes shared points across procs
202964bf7afSMatthew G. Knepley    sfDof - describes shared dofs across procs
203964bf7afSMatthew G. Knepley    sfAdj - describes shared adjacent dofs across procs
204cb1e1211SMatthew G Knepley    ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
20576185916SToby Isaac   (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj
20676185916SToby Isaac        (This is done in DMPlexComputeAnchorAdjacencies())
207cb1e1211SMatthew G Knepley     1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
208cb1e1211SMatthew G Knepley        Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
209cb1e1211SMatthew G Knepley     2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
210cb1e1211SMatthew G Knepley        Create sfAdj connecting rootSectionAdj and leafSectionAdj
211cb1e1211SMatthew G Knepley     3. Visit unowned points on interface, write adjacencies to adj
212cb1e1211SMatthew G Knepley        Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
213cb1e1211SMatthew G Knepley     4. Visit owned points on interface, write adjacencies to rootAdj
214cb1e1211SMatthew G Knepley        Remove redundancy in rootAdj
215cb1e1211SMatthew G Knepley    ** The last two traversals use transitive closure
216cb1e1211SMatthew G Knepley     5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
217cb1e1211SMatthew G Knepley        Allocate memory addressed by sectionAdj (cols)
218cb1e1211SMatthew G Knepley     6. Visit all owned points in the subdomain, insert dof adjacencies into cols
219cb1e1211SMatthew G Knepley    ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
220cb1e1211SMatthew G Knepley   */
221*a9fb6443SMatthew G. Knepley   ierr = DMPlexComputeAnchorAdjacencies(dm, section, sectionGlobal, useCone, useClosure, &anchorSectionAdj, &anchorAdj);CHKERRQ(ierr);
222cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
22376185916SToby Isaac     PetscInt dof, off, d, q, anDof;
22470034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
225cb1e1211SMatthew G Knepley 
226fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
227cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
228cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
229*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
230cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
231cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
232cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
233cb1e1211SMatthew G Knepley 
234cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
235cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
236cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
237cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
238cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
239cb1e1211SMatthew G Knepley       }
240cb1e1211SMatthew G Knepley     }
24176185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
24276185916SToby Isaac     if (anDof) {
24376185916SToby Isaac       for (d = off; d < off+dof; ++d) {
24476185916SToby Isaac         ierr = PetscSectionAddDof(leafSectionAdj, d, anDof);CHKERRQ(ierr);
24576185916SToby Isaac       }
24676185916SToby Isaac     }
247cb1e1211SMatthew G Knepley   }
248cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
249cb1e1211SMatthew G Knepley   if (debug) {
250cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n");CHKERRQ(ierr);
251cb1e1211SMatthew G Knepley     ierr = PetscSectionView(leafSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
252cb1e1211SMatthew G Knepley   }
253cb1e1211SMatthew G Knepley   /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
25447a6104aSMatthew G. Knepley   if (doComm) {
255cb1e1211SMatthew G Knepley     ierr = PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
256cb1e1211SMatthew G Knepley     ierr = PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
257cb1e1211SMatthew G Knepley   }
258cb1e1211SMatthew G Knepley   if (debug) {
259cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n");CHKERRQ(ierr);
260cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
261cb1e1211SMatthew G Knepley   }
262cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
263cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
26476185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof;
265cb1e1211SMatthew G Knepley 
266cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
267cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
268cb1e1211SMatthew G Knepley     if (!dof) continue;
269cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
270cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
271*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
272cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
273cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
274cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
275cb1e1211SMatthew G Knepley 
276cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
277cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
278cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
279cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
280cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
281cb1e1211SMatthew G Knepley       }
282cb1e1211SMatthew G Knepley     }
28376185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
28476185916SToby Isaac     if (anDof) {
28576185916SToby Isaac       for (d = off; d < off+dof; ++d) {
28676185916SToby Isaac         ierr = PetscSectionAddDof(rootSectionAdj, d, anDof);CHKERRQ(ierr);
28776185916SToby Isaac       }
28876185916SToby Isaac     }
289cb1e1211SMatthew G Knepley   }
290cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
291cb1e1211SMatthew G Knepley   if (debug) {
292cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n");CHKERRQ(ierr);
293cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
294cb1e1211SMatthew G Knepley   }
295cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
296cb1e1211SMatthew G Knepley   ierr = PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets);CHKERRQ(ierr);
297cb1e1211SMatthew G Knepley   ierr = PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj);CHKERRQ(ierr);
298cb1e1211SMatthew G Knepley   if (debug) {
299cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency SF for Preallocation:\n");CHKERRQ(ierr);
300cb1e1211SMatthew G Knepley     ierr = PetscSFView(sfAdj, NULL);CHKERRQ(ierr);
301cb1e1211SMatthew G Knepley   }
302cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
303cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
304cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(leafSectionAdj, &adjSize);CHKERRQ(ierr);
3051795a4d1SJed Brown   ierr = PetscCalloc1(adjSize, &adj);CHKERRQ(ierr);
306cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
30776185916SToby Isaac     PetscInt dof, off, d, q, anDof, anOff;
30870034214SMatthew G. Knepley     PetscInt p = leaves[l], numAdj = PETSC_DETERMINE;
309cb1e1211SMatthew G Knepley 
310fb47e2feSMatthew G. Knepley     if ((p < pStart) || (p >= pEnd)) continue;
311cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
312cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
313*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
31476185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
31576185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
316cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
317cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
318cb1e1211SMatthew G Knepley 
319cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(leafSectionAdj, d, &aoff);CHKERRQ(ierr);
320cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
321cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
322cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
323cb1e1211SMatthew G Knepley 
324cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
325cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
326cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
327cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
328cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
329cb1e1211SMatthew G Knepley           adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd;
330cb1e1211SMatthew G Knepley           ++i;
331cb1e1211SMatthew G Knepley         }
332cb1e1211SMatthew G Knepley       }
33376185916SToby Isaac       for (q = 0; q < anDof; q++) {
33476185916SToby Isaac         adj[aoff+i] = anchorAdj[anOff+q];
33576185916SToby Isaac         ++i;
33676185916SToby Isaac       }
337cb1e1211SMatthew G Knepley     }
338cb1e1211SMatthew G Knepley   }
339cb1e1211SMatthew G Knepley   /* Debugging */
340cb1e1211SMatthew G Knepley   if (debug) {
341cb1e1211SMatthew G Knepley     IS tmp;
342cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Leaf adjacency indices\n");CHKERRQ(ierr);
343cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
344cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
3450a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
346cb1e1211SMatthew G Knepley   }
347543482b8SMatthew G. Knepley   /* Gather adjacent indices to root */
348cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(rootSectionAdj, &adjSize);CHKERRQ(ierr);
349785e854fSJed Brown   ierr = PetscMalloc1(adjSize, &rootAdj);CHKERRQ(ierr);
350cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
35147a6104aSMatthew G. Knepley   if (doComm) {
352543482b8SMatthew G. Knepley     const PetscInt *indegree;
353543482b8SMatthew G. Knepley     PetscInt       *remoteadj, radjsize = 0;
354543482b8SMatthew G. Knepley 
355543482b8SMatthew G. Knepley     ierr = PetscSFComputeDegreeBegin(sfAdj, &indegree);CHKERRQ(ierr);
356543482b8SMatthew G. Knepley     ierr = PetscSFComputeDegreeEnd(sfAdj, &indegree);CHKERRQ(ierr);
357543482b8SMatthew G. Knepley     for (p = 0; p < adjSize; ++p) radjsize += indegree[p];
358543482b8SMatthew G. Knepley     ierr = PetscMalloc1(radjsize, &remoteadj);CHKERRQ(ierr);
359543482b8SMatthew G. Knepley     ierr = PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr);
360543482b8SMatthew G. Knepley     ierr = PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr);
3616dba2905SMatthew G. Knepley     for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p-1])) {
362543482b8SMatthew G. Knepley       PetscInt s;
3636dba2905SMatthew G. Knepley       for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l+s] = remoteadj[r];
364543482b8SMatthew G. Knepley     }
365543482b8SMatthew G. Knepley     if (r != radjsize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", r, radjsize);
3666dba2905SMatthew G. Knepley     if (l != adjSize)  SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", l, adjSize);
367543482b8SMatthew G. Knepley     ierr = PetscFree(remoteadj);CHKERRQ(ierr);
368cb1e1211SMatthew G Knepley   }
369cb1e1211SMatthew G Knepley   ierr = PetscSFDestroy(&sfAdj);CHKERRQ(ierr);
370cb1e1211SMatthew G Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
371cb1e1211SMatthew G Knepley   /* Debugging */
372cb1e1211SMatthew G Knepley   if (debug) {
373cb1e1211SMatthew G Knepley     IS tmp;
374cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after gather\n");CHKERRQ(ierr);
375cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
376cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
3770a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
378cb1e1211SMatthew G Knepley   }
379cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
380cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
38176185916SToby Isaac     PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff;
382cb1e1211SMatthew G Knepley 
383cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
384cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
385cb1e1211SMatthew G Knepley     if (!dof) continue;
386cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
387cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
388*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
38976185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
39076185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
391cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
392cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
393cb1e1211SMatthew G Knepley 
394cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
395cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
396cb1e1211SMatthew G Knepley       i    = adof-1;
39776185916SToby Isaac       for (q = 0; q < anDof; q++) {
39876185916SToby Isaac         rootAdj[aoff+i] = anchorAdj[anOff+q];
39976185916SToby Isaac         --i;
40076185916SToby Isaac       }
401cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
402cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
403cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
404cb1e1211SMatthew G Knepley 
405cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
406cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
407cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
408cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
409cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
410cb1e1211SMatthew G Knepley           rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
411cb1e1211SMatthew G Knepley           --i;
412cb1e1211SMatthew G Knepley         }
413cb1e1211SMatthew G Knepley       }
414cb1e1211SMatthew G Knepley     }
415cb1e1211SMatthew G Knepley   }
416cb1e1211SMatthew G Knepley   /* Debugging */
417cb1e1211SMatthew G Knepley   if (debug) {
418cb1e1211SMatthew G Knepley     IS tmp;
419cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices\n");CHKERRQ(ierr);
420cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
421cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
4220a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
423cb1e1211SMatthew G Knepley   }
424cb1e1211SMatthew G Knepley   /* Compress indices */
425cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
426cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
427cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
428cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
429cb1e1211SMatthew G Knepley 
430cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
431cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
432cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
433cb1e1211SMatthew G Knepley     if (!dof) continue;
434cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
435cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
436cb1e1211SMatthew G Knepley     for (d = off; d < off+dof-cdof; ++d) {
437cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
438cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
439cb1e1211SMatthew G Knepley       ierr = PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]);CHKERRQ(ierr);
440cb1e1211SMatthew G Knepley       ierr = PetscSectionSetDof(rootSectionAdj, d, adof);CHKERRQ(ierr);
441cb1e1211SMatthew G Knepley     }
442cb1e1211SMatthew G Knepley   }
443cb1e1211SMatthew G Knepley   /* Debugging */
444cb1e1211SMatthew G Knepley   if (debug) {
445cb1e1211SMatthew G Knepley     IS tmp;
446cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n");CHKERRQ(ierr);
447cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
448cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after compression\n");CHKERRQ(ierr);
449cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
450cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
4510a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
452cb1e1211SMatthew G Knepley   }
453cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
454cb1e1211SMatthew G Knepley   ierr = PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd);CHKERRQ(ierr);
455cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &sectionAdj);CHKERRQ(ierr);
456cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd);CHKERRQ(ierr);
457cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
45876185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof;
459cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
460cb1e1211SMatthew G Knepley 
461cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
462cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
463cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
464cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
465cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
466cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
467cb1e1211SMatthew G Knepley 
468cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
469cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
470cb1e1211SMatthew G Knepley       if (ldof > 0) {
471cb1e1211SMatthew G Knepley         /* We do not own this point */
472cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
473cb1e1211SMatthew G Knepley         ierr = PetscSectionSetDof(sectionAdj, goff+d, rdof);CHKERRQ(ierr);
474cb1e1211SMatthew G Knepley       } else {
475cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
476cb1e1211SMatthew G Knepley       }
477cb1e1211SMatthew G Knepley     }
478cb1e1211SMatthew G Knepley     if (found) continue;
479cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
480cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
481*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
482cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
483cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
484cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof, noff;
485cb1e1211SMatthew G Knepley 
486cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
487cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
488cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
489cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(section, padj, &noff);CHKERRQ(ierr);
490cb1e1211SMatthew G Knepley       for (d = goff; d < goff+dof-cdof; ++d) {
491cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(sectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
492cb1e1211SMatthew G Knepley       }
493cb1e1211SMatthew G Knepley     }
49476185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
49576185916SToby Isaac     if (anDof) {
49676185916SToby Isaac       for (d = goff; d < goff+dof-cdof; ++d) {
49776185916SToby Isaac         ierr = PetscSectionAddDof(sectionAdj, d, anDof);CHKERRQ(ierr);
49876185916SToby Isaac       }
49976185916SToby Isaac     }
500cb1e1211SMatthew G Knepley   }
501cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(sectionAdj);CHKERRQ(ierr);
502cb1e1211SMatthew G Knepley   if (debug) {
503cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation:\n");CHKERRQ(ierr);
504cb1e1211SMatthew G Knepley     ierr = PetscSectionView(sectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
505cb1e1211SMatthew G Knepley   }
506cb1e1211SMatthew G Knepley   /* Get adjacent indices */
507cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(sectionAdj, &numCols);CHKERRQ(ierr);
508785e854fSJed Brown   ierr = PetscMalloc1(numCols, &cols);CHKERRQ(ierr);
509cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
51076185916SToby Isaac     PetscInt  numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff;
511cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
512cb1e1211SMatthew G Knepley 
513cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
514cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
515cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
516cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
517cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
518cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
519cb1e1211SMatthew G Knepley 
520cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
521cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
522cb1e1211SMatthew G Knepley       if (ldof > 0) {
523cb1e1211SMatthew G Knepley         /* We do not own this point */
524cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
525cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
526cb1e1211SMatthew G Knepley 
527cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionAdj, goff+d, &aoff);CHKERRQ(ierr);
528cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(rootSectionAdj, off+d, &roff);CHKERRQ(ierr);
529cb1e1211SMatthew G Knepley         ierr = PetscMemcpy(&cols[aoff], &rootAdj[roff], rdof * sizeof(PetscInt));CHKERRQ(ierr);
530cb1e1211SMatthew G Knepley       } else {
531cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
532cb1e1211SMatthew G Knepley       }
533cb1e1211SMatthew G Knepley     }
534cb1e1211SMatthew G Knepley     if (found) continue;
535*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr);
53676185916SToby Isaac     ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr);
53776185916SToby Isaac     ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr);
538cb1e1211SMatthew G Knepley     for (d = goff; d < goff+dof-cdof; ++d) {
539cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
540cb1e1211SMatthew G Knepley 
541cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(sectionAdj, d, &adof);CHKERRQ(ierr);
542cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(sectionAdj, d, &aoff);CHKERRQ(ierr);
543cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
544cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
545cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
546cb1e1211SMatthew G Knepley         const PetscInt *ncind;
547cb1e1211SMatthew G Knepley 
548cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
549cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
550cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
551cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
552cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintIndices(section, padj, &ncind);CHKERRQ(ierr);
553cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
554cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd, ++i) {
555cb1e1211SMatthew G Knepley           cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
556cb1e1211SMatthew G Knepley         }
557cb1e1211SMatthew G Knepley       }
55876185916SToby Isaac       for (q = 0; q < anDof; q++, i++) {
55976185916SToby Isaac         cols[aoff+i] = anchorAdj[anOff + q];
56076185916SToby Isaac       }
561cb1e1211SMatthew 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);
562cb1e1211SMatthew G Knepley     }
563cb1e1211SMatthew G Knepley   }
56476185916SToby Isaac   ierr = PetscSectionDestroy(&anchorSectionAdj);CHKERRQ(ierr);
565cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&leafSectionAdj);CHKERRQ(ierr);
566cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&rootSectionAdj);CHKERRQ(ierr);
56776185916SToby Isaac   ierr = PetscFree(anchorAdj);CHKERRQ(ierr);
568cb1e1211SMatthew G Knepley   ierr = PetscFree(rootAdj);CHKERRQ(ierr);
56970034214SMatthew G. Knepley   ierr = PetscFree(tmpAdj);CHKERRQ(ierr);
570cb1e1211SMatthew G Knepley   /* Debugging */
571cb1e1211SMatthew G Knepley   if (debug) {
572cb1e1211SMatthew G Knepley     IS tmp;
573cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Column indices\n");CHKERRQ(ierr);
574cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
575cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
5760a7e1f69SMatthew G. Knepley     ierr = ISDestroy(&tmp);CHKERRQ(ierr);
577cb1e1211SMatthew G Knepley   }
578*a9fb6443SMatthew G. Knepley 
579*a9fb6443SMatthew G. Knepley   *sA     = sectionAdj;
580*a9fb6443SMatthew G. Knepley   *colIdx = cols;
581*a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
582*a9fb6443SMatthew G. Knepley }
583*a9fb6443SMatthew G. Knepley 
584*a9fb6443SMatthew G. Knepley #undef __FUNCT__
585*a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexUpdateAllocation_Static"
586*a9fb6443SMatthew G. Knepley PetscErrorCode DMPlexUpdateAllocation_Static(PetscLayout rLayout, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[])
587*a9fb6443SMatthew G. Knepley {
588*a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p;
589*a9fb6443SMatthew G. Knepley   PetscErrorCode ierr;
590*a9fb6443SMatthew G. Knepley 
591*a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
592*a9fb6443SMatthew G. Knepley   /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */
593cb1e1211SMatthew G Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
594*a9fb6443SMatthew 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);
595*a9fb6443SMatthew G. Knepley   if (f >= 0 && bs == 1) {
596*a9fb6443SMatthew G. Knepley     ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
597*a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
598*a9fb6443SMatthew G. Knepley       PetscInt goff, loff, lfoff, fdof, fcdof, rS, rE;
599*a9fb6443SMatthew G. Knepley 
600*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
601*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetOffset(section, p, &loff);CHKERRQ(ierr);
602*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetFieldOffset(section, p, f, &lfoff);CHKERRQ(ierr);
603*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
604*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
605*a9fb6443SMatthew G. Knepley       rS   = goff + (lfoff - loff);
606*a9fb6443SMatthew G. Knepley       rE   = rS + fdof - fcdof;
607*a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
608*a9fb6443SMatthew G. Knepley         PetscInt numCols, cStart, c;
609*a9fb6443SMatthew G. Knepley 
610*a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
611*a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
612*a9fb6443SMatthew G. Knepley         for (c = cStart; c < cStart+numCols; ++c) {
613*a9fb6443SMatthew G. Knepley           if ((cols[c] >= rStart) && (cols[c] < rEnd)) {
614*a9fb6443SMatthew G. Knepley             ++dnz[r-rStart];
615*a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++dnzu[r-rStart];
616*a9fb6443SMatthew G. Knepley           } else {
617*a9fb6443SMatthew G. Knepley             ++onz[r-rStart];
618*a9fb6443SMatthew G. Knepley             if (cols[c] >= r) ++onzu[r-rStart];
619*a9fb6443SMatthew G. Knepley           }
620*a9fb6443SMatthew G. Knepley         }
621*a9fb6443SMatthew G. Knepley       }
622*a9fb6443SMatthew G. Knepley     }
623*a9fb6443SMatthew 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) {
632cb1e1211SMatthew G Knepley         if ((cols[c] >= rStart*bs) && (cols[c] < rEnd*bs)) {
633cb1e1211SMatthew G Knepley           ++dnz[r-rStart];
634cb1e1211SMatthew G Knepley           if (cols[c] >= row) ++dnzu[r-rStart];
635cb1e1211SMatthew G Knepley         } else {
636cb1e1211SMatthew G Knepley           ++onz[r-rStart];
637cb1e1211SMatthew G Knepley           if (cols[c] >= row) ++onzu[r-rStart];
638cb1e1211SMatthew G Knepley         }
639cb1e1211SMatthew G Knepley       }
640cb1e1211SMatthew G Knepley     }
641*a9fb6443SMatthew 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   }
648*a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
649*a9fb6443SMatthew G. Knepley }
650*a9fb6443SMatthew G. Knepley 
651*a9fb6443SMatthew G. Knepley #undef __FUNCT__
652*a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexFillMatrix_Static"
653*a9fb6443SMatthew G. Knepley PetscErrorCode DMPlexFillMatrix_Static(PetscLayout rLayout, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A)
654*a9fb6443SMatthew G. Knepley {
655*a9fb6443SMatthew G. Knepley   PetscScalar   *values;
656*a9fb6443SMatthew G. Knepley   PetscInt       rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0;
657*a9fb6443SMatthew G. Knepley   PetscErrorCode ierr;
658*a9fb6443SMatthew G. Knepley 
659*a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
660*a9fb6443SMatthew G. Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
661*a9fb6443SMatthew G. Knepley   for (r = rStart; r < rEnd; ++r) {
662*a9fb6443SMatthew G. Knepley     ierr      = PetscSectionGetDof(sectionAdj, r, &len);CHKERRQ(ierr);
663*a9fb6443SMatthew G. Knepley     maxRowLen = PetscMax(maxRowLen, len);
664*a9fb6443SMatthew G. Knepley   }
665*a9fb6443SMatthew G. Knepley   ierr = PetscCalloc1(maxRowLen, &values);CHKERRQ(ierr);
666*a9fb6443SMatthew G. Knepley   if (f >=0 && bs == 1) {
667*a9fb6443SMatthew G. Knepley     ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
668*a9fb6443SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
669*a9fb6443SMatthew G. Knepley       PetscInt goff, loff, lfoff, fdof, fcdof, rS, rE;
670*a9fb6443SMatthew G. Knepley 
671*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
672*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetOffset(section, p, &loff);CHKERRQ(ierr);
673*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetFieldOffset(section, p, f, &lfoff);CHKERRQ(ierr);
674*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetFieldDof(section, p, f, &fdof);CHKERRQ(ierr);
675*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr);
676*a9fb6443SMatthew G. Knepley       rS   = goff + (lfoff - loff);
677*a9fb6443SMatthew G. Knepley       rE   = rS + fdof - fcdof;
678*a9fb6443SMatthew G. Knepley       for (r = rS; r < rE; ++r) {
679*a9fb6443SMatthew G. Knepley         PetscInt numCols, cStart, c;
680*a9fb6443SMatthew G. Knepley 
681*a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
682*a9fb6443SMatthew G. Knepley         ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
683*a9fb6443SMatthew G. Knepley         ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr);
684*a9fb6443SMatthew G. Knepley       }
685*a9fb6443SMatthew G. Knepley     }
686*a9fb6443SMatthew G. Knepley   } else {
687*a9fb6443SMatthew G. Knepley     for (r = rStart; r < rEnd; ++r) {
688*a9fb6443SMatthew G. Knepley       PetscInt numCols, cStart;
689*a9fb6443SMatthew G. Knepley 
690*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
691*a9fb6443SMatthew G. Knepley       ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
692*a9fb6443SMatthew G. Knepley       ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr);
693*a9fb6443SMatthew G. Knepley     }
694*a9fb6443SMatthew G. Knepley   }
695*a9fb6443SMatthew G. Knepley   ierr = PetscFree(values);CHKERRQ(ierr);
696*a9fb6443SMatthew G. Knepley   PetscFunctionReturn(0);
697*a9fb6443SMatthew G. Knepley }
698*a9fb6443SMatthew G. Knepley 
699*a9fb6443SMatthew G. Knepley #undef __FUNCT__
700*a9fb6443SMatthew G. Knepley #define __FUNCT__ "DMPlexPreallocateOperator"
701*a9fb6443SMatthew G. Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
702*a9fb6443SMatthew G. Knepley {
703*a9fb6443SMatthew G. Knepley   MPI_Comm       comm;
704*a9fb6443SMatthew G. Knepley   PetscDS        prob;
705*a9fb6443SMatthew G. Knepley   MatType        mtype;
706*a9fb6443SMatthew G. Knepley   PetscSF        sf, sfDof;
707*a9fb6443SMatthew G. Knepley   PetscInt      *remoteOffsets;
708*a9fb6443SMatthew G. Knepley   PetscSection   sectionAdj[4] = {NULL, NULL, NULL, NULL};
709*a9fb6443SMatthew G. Knepley   PetscInt      *cols[4]       = {NULL, NULL, NULL, NULL};
710*a9fb6443SMatthew G. Knepley   PetscBool      useCone, useClosure;
711*a9fb6443SMatthew G. Knepley   PetscInt       Nf, f, idx, locRows, r;
712*a9fb6443SMatthew G. Knepley   PetscLayout    rLayout;
713*a9fb6443SMatthew G. Knepley   PetscBool      isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE;
714*a9fb6443SMatthew G. Knepley   PetscErrorCode ierr;
715*a9fb6443SMatthew G. Knepley 
716*a9fb6443SMatthew G. Knepley   PetscFunctionBegin;
717*a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
718*a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3);
719*a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(sectionGlobal, PETSC_SECTION_CLASSID, 4);
720*a9fb6443SMatthew G. Knepley   PetscValidHeaderSpecific(A, MAT_CLASSID, 9);
721*a9fb6443SMatthew G. Knepley   if (dnz)  PetscValidPointer(dnz,5);  if (onz)  PetscValidPointer(onz,6);
722*a9fb6443SMatthew G. Knepley   if (dnzu) PetscValidPointer(dnzu,7); if (onzu) PetscValidPointer(onzu,8);
723*a9fb6443SMatthew G. Knepley   ierr = DMGetDS(dm, &prob);CHKERRQ(ierr);
724*a9fb6443SMatthew G. Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
725*a9fb6443SMatthew G. Knepley   ierr = PetscOptionsGetBool(NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr);
726*a9fb6443SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
727*a9fb6443SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr);
728*a9fb6443SMatthew G. Knepley   /* Create dof SF based on point SF */
729*a9fb6443SMatthew G. Knepley   if (debug) {
730*a9fb6443SMatthew G. Knepley     PetscSF sf;
731*a9fb6443SMatthew G. Knepley 
732*a9fb6443SMatthew G. Knepley     ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
733*a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Input Section for Preallocation:\n");CHKERRQ(ierr);
734*a9fb6443SMatthew G. Knepley     ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
735*a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Input Global Section for Preallocation:\n");CHKERRQ(ierr);
736*a9fb6443SMatthew G. Knepley     ierr = PetscSectionView(sectionGlobal, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
737*a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Input SF for Preallocation:\n");CHKERRQ(ierr);
738*a9fb6443SMatthew G. Knepley     ierr = PetscSFView(sf, NULL);CHKERRQ(ierr);
739*a9fb6443SMatthew G. Knepley   }
740*a9fb6443SMatthew G. Knepley   ierr = PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets);CHKERRQ(ierr);
741*a9fb6443SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof);CHKERRQ(ierr);
742*a9fb6443SMatthew G. Knepley   if (debug) {
743*a9fb6443SMatthew G. Knepley     ierr = PetscPrintf(comm, "Dof SF for Preallocation:\n");CHKERRQ(ierr);
744*a9fb6443SMatthew G. Knepley     ierr = PetscSFView(sfDof, NULL);CHKERRQ(ierr);
745*a9fb6443SMatthew G. Knepley   }
746*a9fb6443SMatthew G. Knepley   /* Create allocation vectors from adjacency graph */
747*a9fb6443SMatthew G. Knepley   ierr = MatGetLocalSize(A, &locRows, NULL);CHKERRQ(ierr);
748*a9fb6443SMatthew G. Knepley   ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr);
749*a9fb6443SMatthew G. Knepley   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
750*a9fb6443SMatthew G. Knepley   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
751*a9fb6443SMatthew G. Knepley   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
752*a9fb6443SMatthew G. Knepley   /* There are 4 types of adjacency */
753*a9fb6443SMatthew G. Knepley   ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr);
754*a9fb6443SMatthew G. Knepley   if (Nf < 1 || bs > 1) {
755*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
756*a9fb6443SMatthew G. Knepley     ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
757*a9fb6443SMatthew G. Knepley     idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
758*a9fb6443SMatthew G. Knepley     ierr = DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, section, sectionGlobal, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]);CHKERRQ(ierr);
759*a9fb6443SMatthew G. Knepley     ierr = DMPlexUpdateAllocation_Static(rLayout, bs, section, sectionGlobal, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu);CHKERRQ(ierr);
760*a9fb6443SMatthew G. Knepley   } else {
761*a9fb6443SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
762*a9fb6443SMatthew G. Knepley       ierr = PetscDSGetAdjacency(prob, f, &useCone, &useClosure);CHKERRQ(ierr);
763*a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
764*a9fb6443SMatthew G. Knepley       if (!sectionAdj[idx]) {ierr = DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, section, sectionGlobal, useCone, useClosure, PETSC_TRUE, &sectionAdj[idx], &cols[idx]);CHKERRQ(ierr);}
765*a9fb6443SMatthew G. Knepley       ierr = DMPlexUpdateAllocation_Static(rLayout, bs, section, sectionGlobal, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu);CHKERRQ(ierr);
766*a9fb6443SMatthew G. Knepley     }
767*a9fb6443SMatthew G. Knepley   }
768*a9fb6443SMatthew 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) {
780*a9fb6443SMatthew G. Knepley     if (Nf < 1 || bs > 1) {
781*a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
782*a9fb6443SMatthew G. Knepley       ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
783*a9fb6443SMatthew G. Knepley       idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
784*a9fb6443SMatthew G. Knepley       ierr = DMPlexFillMatrix_Static(rLayout, bs, section, sectionGlobal, -1, sectionAdj[idx], cols[idx], A);CHKERRQ(ierr);
785*a9fb6443SMatthew G. Knepley     } else {
786*a9fb6443SMatthew G. Knepley       for (f = 0; f < Nf; ++f) {
787*a9fb6443SMatthew G. Knepley         ierr = PetscDSGetAdjacency(prob, f, &useCone, &useClosure);CHKERRQ(ierr);
788*a9fb6443SMatthew G. Knepley         idx  = (useCone ? 1 : 0) + (useClosure ? 2 : 0);
789*a9fb6443SMatthew G. Knepley         ierr = DMPlexFillMatrix_Static(rLayout, bs, section, sectionGlobal, 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   }
795*a9fb6443SMatthew G. Knepley   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
796*a9fb6443SMatthew 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