xref: /petsc/src/dm/impls/plex/plexpreallocate.c (revision cb1e1211e79a33490225712e6f7f526b093c4371)
1*cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*cb1e1211SMatthew G Knepley #include <petsc-private/isimpl.h>
3*cb1e1211SMatthew G Knepley 
4*cb1e1211SMatthew G Knepley #undef __FUNCT__
5*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal"
6*cb1e1211SMatthew G Knepley static PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useClosure, const PetscInt *tmpClosure, PetscInt *adjSize, PetscInt adj[])
7*cb1e1211SMatthew G Knepley {
8*cb1e1211SMatthew G Knepley   const PetscInt *star  = tmpClosure;
9*cb1e1211SMatthew G Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, starSize, s;
10*cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
11*cb1e1211SMatthew G Knepley 
12*cb1e1211SMatthew G Knepley   PetscFunctionBegin;
13*cb1e1211SMatthew G Knepley   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, (PetscInt**) &star);CHKERRQ(ierr);
14*cb1e1211SMatthew G Knepley   for (s = 2; s < starSize*2; s += 2) {
15*cb1e1211SMatthew G Knepley     const PetscInt *closure = NULL;
16*cb1e1211SMatthew G Knepley     PetscInt        closureSize, c, q;
17*cb1e1211SMatthew G Knepley 
18*cb1e1211SMatthew G Knepley     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
19*cb1e1211SMatthew G Knepley     for (c = 0; c < closureSize*2; c += 2) {
20*cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
21*cb1e1211SMatthew G Knepley         if (closure[c] == adj[q]) break;
22*cb1e1211SMatthew G Knepley       }
23*cb1e1211SMatthew G Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
24*cb1e1211SMatthew G Knepley     }
25*cb1e1211SMatthew G Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
26*cb1e1211SMatthew G Knepley   }
27*cb1e1211SMatthew G Knepley   *adjSize = numAdj;
28*cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
29*cb1e1211SMatthew G Knepley }
30*cb1e1211SMatthew G Knepley 
31*cb1e1211SMatthew G Knepley #undef __FUNCT__
32*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetPreallocationCenterDimension"
33*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetPreallocationCenterDimension(DM dm, PetscInt preallocCenterDim)
34*cb1e1211SMatthew G Knepley {
35*cb1e1211SMatthew G Knepley   DM_Plex *mesh = (DM_Plex*) dm->data;
36*cb1e1211SMatthew G Knepley 
37*cb1e1211SMatthew G Knepley   PetscFunctionBegin;
38*cb1e1211SMatthew G Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
39*cb1e1211SMatthew G Knepley   mesh->preallocCenterDim = preallocCenterDim;
40*cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
41*cb1e1211SMatthew G Knepley }
42*cb1e1211SMatthew G Knepley 
43*cb1e1211SMatthew G Knepley #undef __FUNCT__
44*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator"
45*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix)
46*cb1e1211SMatthew G Knepley {
47*cb1e1211SMatthew G Knepley   DM_Plex           *mesh = (DM_Plex*) dm->data;
48*cb1e1211SMatthew G Knepley   MPI_Comm           comm;
49*cb1e1211SMatthew G Knepley   PetscSF            sf, sfDof, sfAdj;
50*cb1e1211SMatthew G Knepley   PetscSection       leafSectionAdj, rootSectionAdj, sectionAdj;
51*cb1e1211SMatthew G Knepley   PetscInt           nleaves, l, p;
52*cb1e1211SMatthew G Knepley   const PetscInt    *leaves;
53*cb1e1211SMatthew G Knepley   const PetscSFNode *remotes;
54*cb1e1211SMatthew G Knepley   PetscInt           dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols;
55*cb1e1211SMatthew G Knepley   PetscInt          *tmpClosure, *tmpAdj, *adj, *rootAdj, *cols, *remoteOffsets;
56*cb1e1211SMatthew G Knepley   PetscInt           depth, maxConeSize, maxSupportSize, maxClosureSize, maxAdjSize, adjSize;
57*cb1e1211SMatthew G Knepley   PetscLayout        rLayout;
58*cb1e1211SMatthew G Knepley   PetscInt           locRows, rStart, rEnd, r;
59*cb1e1211SMatthew G Knepley   PetscMPIInt        size;
60*cb1e1211SMatthew G Knepley   PetscBool          useClosure, debug = PETSC_FALSE;
61*cb1e1211SMatthew G Knepley   PetscErrorCode     ierr;
62*cb1e1211SMatthew G Knepley 
63*cb1e1211SMatthew G Knepley   PetscFunctionBegin;
64*cb1e1211SMatthew G Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
65*cb1e1211SMatthew G Knepley   ierr = PetscOptionsGetBool(NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr);
66*cb1e1211SMatthew G Knepley   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
67*cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
68*cb1e1211SMatthew G Knepley   ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr);
69*cb1e1211SMatthew G Knepley   /* Create dof SF based on point SF */
70*cb1e1211SMatthew G Knepley   if (debug) {
71*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Input Section for Preallocation:\n");CHKERRQ(ierr);
72*cb1e1211SMatthew G Knepley     ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
73*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Input Global Section for Preallocation:\n");CHKERRQ(ierr);
74*cb1e1211SMatthew G Knepley     ierr = PetscSectionView(sectionGlobal, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
75*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Input SF for Preallocation:\n");CHKERRQ(ierr);
76*cb1e1211SMatthew G Knepley     ierr = PetscSFView(sf, NULL);CHKERRQ(ierr);
77*cb1e1211SMatthew G Knepley   }
78*cb1e1211SMatthew G Knepley   ierr = PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets);CHKERRQ(ierr);
79*cb1e1211SMatthew G Knepley   ierr = PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof);CHKERRQ(ierr);
80*cb1e1211SMatthew G Knepley   if (debug) {
81*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Dof SF for Preallocation:\n");CHKERRQ(ierr);
82*cb1e1211SMatthew G Knepley     ierr = PetscSFView(sfDof, NULL);CHKERRQ(ierr);
83*cb1e1211SMatthew G Knepley   }
84*cb1e1211SMatthew G Knepley   /* Create section for dof adjacency (dof ==> # adj dof) */
85*cb1e1211SMatthew G Knepley   /*   FEM:   Two points p and q are adjacent if q \in closure(star(p)), preallocCenterDim = dim   */
86*cb1e1211SMatthew G Knepley   /*   FVM:   Two points p and q are adjacent if q \in star(cone(p)),    preallocCenterDim = dim-1 */
87*cb1e1211SMatthew G Knepley   /*   FVM++: Two points p and q are adjacent if q \in star(closure(p)), preallocCenterDim = 0     */
88*cb1e1211SMatthew G Knepley   if (mesh->preallocCenterDim == dim) {
89*cb1e1211SMatthew G Knepley     useClosure = PETSC_FALSE;
90*cb1e1211SMatthew G Knepley   } else if (mesh->preallocCenterDim == 0) {
91*cb1e1211SMatthew G Knepley     useClosure = PETSC_TRUE;
92*cb1e1211SMatthew G Knepley   } else SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Do not support preallocation with center points of dimension %d", mesh->preallocCenterDim);
93*cb1e1211SMatthew G Knepley 
94*cb1e1211SMatthew G Knepley   ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
95*cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(section, &numDof);CHKERRQ(ierr);
96*cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &leafSectionAdj);CHKERRQ(ierr);
97*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(leafSectionAdj, 0, numDof);CHKERRQ(ierr);
98*cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &rootSectionAdj);CHKERRQ(ierr);
99*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(rootSectionAdj, 0, numDof);CHKERRQ(ierr);
100*cb1e1211SMatthew G Knepley   /*   Fill in the ghost dofs on the interface */
101*cb1e1211SMatthew G Knepley   ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes);CHKERRQ(ierr);
102*cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
103*cb1e1211SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
104*cb1e1211SMatthew G Knepley 
105*cb1e1211SMatthew G Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth),PetscPowInt(mesh->maxSupportSize,depth)) + 2;
106*cb1e1211SMatthew G Knepley   maxAdjSize     = PetscPowInt(mesh->maxConeSize,depth) * PetscPowInt(mesh->maxSupportSize,depth) + 1;
107*cb1e1211SMatthew G Knepley 
108*cb1e1211SMatthew G Knepley   ierr = PetscMalloc2(maxClosureSize,PetscInt,&tmpClosure,maxAdjSize,PetscInt,&tmpAdj);CHKERRQ(ierr);
109*cb1e1211SMatthew G Knepley 
110*cb1e1211SMatthew G Knepley   /*
111*cb1e1211SMatthew G Knepley    ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point.
112*cb1e1211SMatthew G Knepley     1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj
113*cb1e1211SMatthew G Knepley        Reduce those counts to rootSectionAdj (now redundantly counting some interface points)
114*cb1e1211SMatthew G Knepley     2. Visit owned points on interface, count adjacencies placing in rootSectionAdj
115*cb1e1211SMatthew G Knepley        Create sfAdj connecting rootSectionAdj and leafSectionAdj
116*cb1e1211SMatthew G Knepley     3. Visit unowned points on interface, write adjacencies to adj
117*cb1e1211SMatthew G Knepley        Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies)
118*cb1e1211SMatthew G Knepley     4. Visit owned points on interface, write adjacencies to rootAdj
119*cb1e1211SMatthew G Knepley        Remove redundancy in rootAdj
120*cb1e1211SMatthew G Knepley    ** The last two traversals use transitive closure
121*cb1e1211SMatthew G Knepley     5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj)
122*cb1e1211SMatthew G Knepley        Allocate memory addressed by sectionAdj (cols)
123*cb1e1211SMatthew G Knepley     6. Visit all owned points in the subdomain, insert dof adjacencies into cols
124*cb1e1211SMatthew G Knepley    ** Knowing all the column adjacencies, check ownership and sum into dnz and onz
125*cb1e1211SMatthew G Knepley   */
126*cb1e1211SMatthew G Knepley 
127*cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
128*cb1e1211SMatthew G Knepley     PetscInt dof, off, d, q;
129*cb1e1211SMatthew G Knepley     PetscInt p = leaves[l], numAdj = maxAdjSize;
130*cb1e1211SMatthew G Knepley 
131*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
132*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
133*cb1e1211SMatthew G Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useClosure, tmpClosure, &numAdj, tmpAdj);CHKERRQ(ierr);
134*cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
135*cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
136*cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
137*cb1e1211SMatthew G Knepley 
138*cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
139*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
140*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
141*cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
142*cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
143*cb1e1211SMatthew G Knepley       }
144*cb1e1211SMatthew G Knepley     }
145*cb1e1211SMatthew G Knepley   }
146*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
147*cb1e1211SMatthew G Knepley   if (debug) {
148*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n");CHKERRQ(ierr);
149*cb1e1211SMatthew G Knepley     ierr = PetscSectionView(leafSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
150*cb1e1211SMatthew G Knepley   }
151*cb1e1211SMatthew G Knepley   /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */
152*cb1e1211SMatthew G Knepley   if (size > 1) {
153*cb1e1211SMatthew G Knepley     ierr = PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
154*cb1e1211SMatthew G Knepley     ierr = PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr);
155*cb1e1211SMatthew G Knepley   }
156*cb1e1211SMatthew G Knepley   if (debug) {
157*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n");CHKERRQ(ierr);
158*cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
159*cb1e1211SMatthew G Knepley   }
160*cb1e1211SMatthew G Knepley   /* Add in local adjacency sizes for owned dofs on interface (roots) */
161*cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
162*cb1e1211SMatthew G Knepley     PetscInt numAdj = maxAdjSize, adof, dof, off, d, q;
163*cb1e1211SMatthew G Knepley 
164*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
165*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
166*cb1e1211SMatthew G Knepley     if (!dof) continue;
167*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
168*cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
169*cb1e1211SMatthew G Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useClosure, tmpClosure, &numAdj, tmpAdj);CHKERRQ(ierr);
170*cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
171*cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
172*cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof;
173*cb1e1211SMatthew G Knepley 
174*cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
175*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
176*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
177*cb1e1211SMatthew G Knepley       for (d = off; d < off+dof; ++d) {
178*cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
179*cb1e1211SMatthew G Knepley       }
180*cb1e1211SMatthew G Knepley     }
181*cb1e1211SMatthew G Knepley   }
182*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
183*cb1e1211SMatthew G Knepley   if (debug) {
184*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n");CHKERRQ(ierr);
185*cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
186*cb1e1211SMatthew G Knepley   }
187*cb1e1211SMatthew G Knepley   /* Create adj SF based on dof SF */
188*cb1e1211SMatthew G Knepley   ierr = PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets);CHKERRQ(ierr);
189*cb1e1211SMatthew G Knepley   ierr = PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj);CHKERRQ(ierr);
190*cb1e1211SMatthew G Knepley   if (debug) {
191*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency SF for Preallocation:\n");CHKERRQ(ierr);
192*cb1e1211SMatthew G Knepley     ierr = PetscSFView(sfAdj, NULL);CHKERRQ(ierr);
193*cb1e1211SMatthew G Knepley   }
194*cb1e1211SMatthew G Knepley   ierr = PetscSFDestroy(&sfDof);CHKERRQ(ierr);
195*cb1e1211SMatthew G Knepley   /* Create leaf adjacency */
196*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr);
197*cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(leafSectionAdj, &adjSize);CHKERRQ(ierr);
198*cb1e1211SMatthew G Knepley   ierr = PetscMalloc(adjSize * sizeof(PetscInt), &adj);CHKERRQ(ierr);
199*cb1e1211SMatthew G Knepley   ierr = PetscMemzero(adj, adjSize * sizeof(PetscInt));CHKERRQ(ierr);
200*cb1e1211SMatthew G Knepley   for (l = 0; l < nleaves; ++l) {
201*cb1e1211SMatthew G Knepley     PetscInt dof, off, d, q;
202*cb1e1211SMatthew G Knepley     PetscInt p = leaves[l], numAdj = maxAdjSize;
203*cb1e1211SMatthew G Knepley 
204*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
205*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
206*cb1e1211SMatthew G Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useClosure, tmpClosure, &numAdj, tmpAdj);CHKERRQ(ierr);
207*cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
208*cb1e1211SMatthew G Knepley       PetscInt aoff, i = 0;
209*cb1e1211SMatthew G Knepley 
210*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(leafSectionAdj, d, &aoff);CHKERRQ(ierr);
211*cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
212*cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
213*cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
214*cb1e1211SMatthew G Knepley 
215*cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
216*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
217*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
218*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
219*cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
220*cb1e1211SMatthew G Knepley           adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd;
221*cb1e1211SMatthew G Knepley           ++i;
222*cb1e1211SMatthew G Knepley         }
223*cb1e1211SMatthew G Knepley       }
224*cb1e1211SMatthew G Knepley     }
225*cb1e1211SMatthew G Knepley   }
226*cb1e1211SMatthew G Knepley   /* Debugging */
227*cb1e1211SMatthew G Knepley   if (debug) {
228*cb1e1211SMatthew G Knepley     IS tmp;
229*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Leaf adjacency indices\n");CHKERRQ(ierr);
230*cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
231*cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
232*cb1e1211SMatthew G Knepley   }
233*cb1e1211SMatthew G Knepley   /* Gather adjacenct indices to root */
234*cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(rootSectionAdj, &adjSize);CHKERRQ(ierr);
235*cb1e1211SMatthew G Knepley   ierr = PetscMalloc(adjSize * sizeof(PetscInt), &rootAdj);CHKERRQ(ierr);
236*cb1e1211SMatthew G Knepley   for (r = 0; r < adjSize; ++r) rootAdj[r] = -1;
237*cb1e1211SMatthew G Knepley   if (size > 1) {
238*cb1e1211SMatthew G Knepley     ierr = PetscSFGatherBegin(sfAdj, MPIU_INT, adj, rootAdj);CHKERRQ(ierr);
239*cb1e1211SMatthew G Knepley     ierr = PetscSFGatherEnd(sfAdj, MPIU_INT, adj, rootAdj);CHKERRQ(ierr);
240*cb1e1211SMatthew G Knepley   }
241*cb1e1211SMatthew G Knepley   ierr = PetscSFDestroy(&sfAdj);CHKERRQ(ierr);
242*cb1e1211SMatthew G Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
243*cb1e1211SMatthew G Knepley   /* Debugging */
244*cb1e1211SMatthew G Knepley   if (debug) {
245*cb1e1211SMatthew G Knepley     IS tmp;
246*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after gather\n");CHKERRQ(ierr);
247*cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
248*cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
249*cb1e1211SMatthew G Knepley   }
250*cb1e1211SMatthew G Knepley   /* Add in local adjacency indices for owned dofs on interface (roots) */
251*cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
252*cb1e1211SMatthew G Knepley     PetscInt numAdj = maxAdjSize, adof, dof, off, d, q;
253*cb1e1211SMatthew G Knepley 
254*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
255*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
256*cb1e1211SMatthew G Knepley     if (!dof) continue;
257*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
258*cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
259*cb1e1211SMatthew G Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useClosure, tmpClosure, &numAdj, tmpAdj);CHKERRQ(ierr);
260*cb1e1211SMatthew G Knepley     for (d = off; d < off+dof; ++d) {
261*cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i;
262*cb1e1211SMatthew G Knepley 
263*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
264*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
265*cb1e1211SMatthew G Knepley       i    = adof-1;
266*cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
267*cb1e1211SMatthew G Knepley         const PetscInt padj = tmpAdj[q];
268*cb1e1211SMatthew G Knepley         PetscInt ndof, ncdof, ngoff, nd;
269*cb1e1211SMatthew G Knepley 
270*cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
271*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
272*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
273*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
274*cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd) {
275*cb1e1211SMatthew G Knepley           rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
276*cb1e1211SMatthew G Knepley           --i;
277*cb1e1211SMatthew G Knepley         }
278*cb1e1211SMatthew G Knepley       }
279*cb1e1211SMatthew G Knepley     }
280*cb1e1211SMatthew G Knepley   }
281*cb1e1211SMatthew G Knepley   /* Debugging */
282*cb1e1211SMatthew G Knepley   if (debug) {
283*cb1e1211SMatthew G Knepley     IS tmp;
284*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices\n");CHKERRQ(ierr);
285*cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
286*cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
287*cb1e1211SMatthew G Knepley   }
288*cb1e1211SMatthew G Knepley   /* Compress indices */
289*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr);
290*cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
291*cb1e1211SMatthew G Knepley     PetscInt dof, cdof, off, d;
292*cb1e1211SMatthew G Knepley     PetscInt adof, aoff;
293*cb1e1211SMatthew G Knepley 
294*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
295*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
296*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
297*cb1e1211SMatthew G Knepley     if (!dof) continue;
298*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr);
299*cb1e1211SMatthew G Knepley     if (adof <= 0) continue;
300*cb1e1211SMatthew G Knepley     for (d = off; d < off+dof-cdof; ++d) {
301*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr);
302*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr);
303*cb1e1211SMatthew G Knepley       ierr = PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]);CHKERRQ(ierr);
304*cb1e1211SMatthew G Knepley       ierr = PetscSectionSetDof(rootSectionAdj, d, adof);CHKERRQ(ierr);
305*cb1e1211SMatthew G Knepley     }
306*cb1e1211SMatthew G Knepley   }
307*cb1e1211SMatthew G Knepley   /* Debugging */
308*cb1e1211SMatthew G Knepley   if (debug) {
309*cb1e1211SMatthew G Knepley     IS tmp;
310*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n");CHKERRQ(ierr);
311*cb1e1211SMatthew G Knepley     ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
312*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Root adjacency indices after compression\n");CHKERRQ(ierr);
313*cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
314*cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
315*cb1e1211SMatthew G Knepley   }
316*cb1e1211SMatthew G Knepley   /* Build adjacency section: Maps global indices to sets of adjacent global indices */
317*cb1e1211SMatthew G Knepley   ierr = PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd);CHKERRQ(ierr);
318*cb1e1211SMatthew G Knepley   ierr = PetscSectionCreate(comm, &sectionAdj);CHKERRQ(ierr);
319*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd);CHKERRQ(ierr);
320*cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
321*cb1e1211SMatthew G Knepley     PetscInt  numAdj = maxAdjSize, dof, cdof, off, goff, d, q;
322*cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
323*cb1e1211SMatthew G Knepley 
324*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
325*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
326*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
327*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
328*cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
329*cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
330*cb1e1211SMatthew G Knepley 
331*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
332*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
333*cb1e1211SMatthew G Knepley       if (ldof > 0) {
334*cb1e1211SMatthew G Knepley         /* We do not own this point */
335*cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
336*cb1e1211SMatthew G Knepley         ierr = PetscSectionSetDof(sectionAdj, goff+d, rdof);CHKERRQ(ierr);
337*cb1e1211SMatthew G Knepley       } else {
338*cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
339*cb1e1211SMatthew G Knepley       }
340*cb1e1211SMatthew G Knepley     }
341*cb1e1211SMatthew G Knepley     if (found) continue;
342*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
343*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
344*cb1e1211SMatthew G Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useClosure, tmpClosure, &numAdj, tmpAdj);CHKERRQ(ierr);
345*cb1e1211SMatthew G Knepley     for (q = 0; q < numAdj; ++q) {
346*cb1e1211SMatthew G Knepley       const PetscInt padj = tmpAdj[q];
347*cb1e1211SMatthew G Knepley       PetscInt ndof, ncdof, noff;
348*cb1e1211SMatthew G Knepley 
349*cb1e1211SMatthew G Knepley       if ((padj < pStart) || (padj >= pEnd)) continue;
350*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
351*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
352*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(section, padj, &noff);CHKERRQ(ierr);
353*cb1e1211SMatthew G Knepley       for (d = goff; d < goff+dof-cdof; ++d) {
354*cb1e1211SMatthew G Knepley         ierr = PetscSectionAddDof(sectionAdj, d, ndof-ncdof);CHKERRQ(ierr);
355*cb1e1211SMatthew G Knepley       }
356*cb1e1211SMatthew G Knepley     }
357*cb1e1211SMatthew G Knepley   }
358*cb1e1211SMatthew G Knepley   ierr = PetscSectionSetUp(sectionAdj);CHKERRQ(ierr);
359*cb1e1211SMatthew G Knepley   if (debug) {
360*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Adjacency Section for Preallocation:\n");CHKERRQ(ierr);
361*cb1e1211SMatthew G Knepley     ierr = PetscSectionView(sectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
362*cb1e1211SMatthew G Knepley   }
363*cb1e1211SMatthew G Knepley   /* Get adjacent indices */
364*cb1e1211SMatthew G Knepley   ierr = PetscSectionGetStorageSize(sectionAdj, &numCols);CHKERRQ(ierr);
365*cb1e1211SMatthew G Knepley   ierr = PetscMalloc(numCols * sizeof(PetscInt), &cols);CHKERRQ(ierr);
366*cb1e1211SMatthew G Knepley   for (p = pStart; p < pEnd; ++p) {
367*cb1e1211SMatthew G Knepley     PetscInt  numAdj = maxAdjSize, dof, cdof, off, goff, d, q;
368*cb1e1211SMatthew G Knepley     PetscBool found  = PETSC_TRUE;
369*cb1e1211SMatthew G Knepley 
370*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
371*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr);
372*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
373*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr);
374*cb1e1211SMatthew G Knepley     for (d = 0; d < dof-cdof; ++d) {
375*cb1e1211SMatthew G Knepley       PetscInt ldof, rdof;
376*cb1e1211SMatthew G Knepley 
377*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr);
378*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr);
379*cb1e1211SMatthew G Knepley       if (ldof > 0) {
380*cb1e1211SMatthew G Knepley         /* We do not own this point */
381*cb1e1211SMatthew G Knepley       } else if (rdof > 0) {
382*cb1e1211SMatthew G Knepley         PetscInt aoff, roff;
383*cb1e1211SMatthew G Knepley 
384*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionAdj, goff+d, &aoff);CHKERRQ(ierr);
385*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(rootSectionAdj, off+d, &roff);CHKERRQ(ierr);
386*cb1e1211SMatthew G Knepley         ierr = PetscMemcpy(&cols[aoff], &rootAdj[roff], rdof * sizeof(PetscInt));CHKERRQ(ierr);
387*cb1e1211SMatthew G Knepley       } else {
388*cb1e1211SMatthew G Knepley         found = PETSC_FALSE;
389*cb1e1211SMatthew G Knepley       }
390*cb1e1211SMatthew G Knepley     }
391*cb1e1211SMatthew G Knepley     if (found) continue;
392*cb1e1211SMatthew G Knepley     ierr = DMPlexGetAdjacency_Internal(dm, p, useClosure, tmpClosure, &numAdj, tmpAdj);CHKERRQ(ierr);
393*cb1e1211SMatthew G Knepley     for (d = goff; d < goff+dof-cdof; ++d) {
394*cb1e1211SMatthew G Knepley       PetscInt adof, aoff, i = 0;
395*cb1e1211SMatthew G Knepley 
396*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(sectionAdj, d, &adof);CHKERRQ(ierr);
397*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(sectionAdj, d, &aoff);CHKERRQ(ierr);
398*cb1e1211SMatthew G Knepley       for (q = 0; q < numAdj; ++q) {
399*cb1e1211SMatthew G Knepley         const PetscInt  padj = tmpAdj[q];
400*cb1e1211SMatthew G Knepley         PetscInt        ndof, ncdof, ngoff, nd;
401*cb1e1211SMatthew G Knepley         const PetscInt *ncind;
402*cb1e1211SMatthew G Knepley 
403*cb1e1211SMatthew G Knepley         /* Adjacent points may not be in the section chart */
404*cb1e1211SMatthew G Knepley         if ((padj < pStart) || (padj >= pEnd)) continue;
405*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr);
406*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr);
407*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetConstraintIndices(section, padj, &ncind);CHKERRQ(ierr);
408*cb1e1211SMatthew G Knepley         ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr);
409*cb1e1211SMatthew G Knepley         for (nd = 0; nd < ndof-ncdof; ++nd, ++i) {
410*cb1e1211SMatthew G Knepley           cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd;
411*cb1e1211SMatthew G Knepley         }
412*cb1e1211SMatthew G Knepley       }
413*cb1e1211SMatthew 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);
414*cb1e1211SMatthew G Knepley     }
415*cb1e1211SMatthew G Knepley   }
416*cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&leafSectionAdj);CHKERRQ(ierr);
417*cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&rootSectionAdj);CHKERRQ(ierr);
418*cb1e1211SMatthew G Knepley   ierr = PetscFree(rootAdj);CHKERRQ(ierr);
419*cb1e1211SMatthew G Knepley   ierr = PetscFree2(tmpClosure, tmpAdj);CHKERRQ(ierr);
420*cb1e1211SMatthew G Knepley   /* Debugging */
421*cb1e1211SMatthew G Knepley   if (debug) {
422*cb1e1211SMatthew G Knepley     IS tmp;
423*cb1e1211SMatthew G Knepley     ierr = PetscPrintf(comm, "Column indices\n");CHKERRQ(ierr);
424*cb1e1211SMatthew G Knepley     ierr = ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr);
425*cb1e1211SMatthew G Knepley     ierr = ISView(tmp, NULL);CHKERRQ(ierr);
426*cb1e1211SMatthew G Knepley   }
427*cb1e1211SMatthew G Knepley   /* Create allocation vectors from adjacency graph */
428*cb1e1211SMatthew G Knepley   ierr = MatGetLocalSize(A, &locRows, NULL);CHKERRQ(ierr);
429*cb1e1211SMatthew G Knepley   ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)A), &rLayout);CHKERRQ(ierr);
430*cb1e1211SMatthew G Knepley   ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr);
431*cb1e1211SMatthew G Knepley   ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr);
432*cb1e1211SMatthew G Knepley   ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr);
433*cb1e1211SMatthew G Knepley   ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr);
434*cb1e1211SMatthew G Knepley   ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr);
435*cb1e1211SMatthew G Knepley   /* Only loop over blocks of rows */
436*cb1e1211SMatthew G Knepley   if (rStart%bs || rEnd%bs) SETERRQ3(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Invalid layout [%d, %d) for matrix, must be divisible by block size %d", rStart, rEnd, bs);
437*cb1e1211SMatthew G Knepley   for (r = rStart/bs; r < rEnd/bs; ++r) {
438*cb1e1211SMatthew G Knepley     const PetscInt row = r*bs;
439*cb1e1211SMatthew G Knepley     PetscInt       numCols, cStart, c;
440*cb1e1211SMatthew G Knepley 
441*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetDof(sectionAdj, row, &numCols);CHKERRQ(ierr);
442*cb1e1211SMatthew G Knepley     ierr = PetscSectionGetOffset(sectionAdj, row, &cStart);CHKERRQ(ierr);
443*cb1e1211SMatthew G Knepley     for (c = cStart; c < cStart+numCols; ++c) {
444*cb1e1211SMatthew G Knepley       if ((cols[c] >= rStart*bs) && (cols[c] < rEnd*bs)) {
445*cb1e1211SMatthew G Knepley         ++dnz[r-rStart];
446*cb1e1211SMatthew G Knepley         if (cols[c] >= row) ++dnzu[r-rStart];
447*cb1e1211SMatthew G Knepley       } else {
448*cb1e1211SMatthew G Knepley         ++onz[r-rStart];
449*cb1e1211SMatthew G Knepley         if (cols[c] >= row) ++onzu[r-rStart];
450*cb1e1211SMatthew G Knepley       }
451*cb1e1211SMatthew G Knepley     }
452*cb1e1211SMatthew G Knepley   }
453*cb1e1211SMatthew G Knepley   if (bs > 1) {
454*cb1e1211SMatthew G Knepley     for (r = 0; r < locRows/bs; ++r) {
455*cb1e1211SMatthew G Knepley       dnz[r]  /= bs;
456*cb1e1211SMatthew G Knepley       onz[r]  /= bs;
457*cb1e1211SMatthew G Knepley       dnzu[r] /= bs;
458*cb1e1211SMatthew G Knepley       onzu[r] /= bs;
459*cb1e1211SMatthew G Knepley     }
460*cb1e1211SMatthew G Knepley   }
461*cb1e1211SMatthew G Knepley   /* Set matrix pattern */
462*cb1e1211SMatthew G Knepley   ierr = MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu);CHKERRQ(ierr);
463*cb1e1211SMatthew G Knepley   ierr = MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
464*cb1e1211SMatthew G Knepley   /* Fill matrix with zeros */
465*cb1e1211SMatthew G Knepley   if (fillMatrix) {
466*cb1e1211SMatthew G Knepley     PetscScalar *values;
467*cb1e1211SMatthew G Knepley     PetscInt     maxRowLen = 0;
468*cb1e1211SMatthew G Knepley 
469*cb1e1211SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
470*cb1e1211SMatthew G Knepley       PetscInt len;
471*cb1e1211SMatthew G Knepley 
472*cb1e1211SMatthew G Knepley       ierr      = PetscSectionGetDof(sectionAdj, r, &len);CHKERRQ(ierr);
473*cb1e1211SMatthew G Knepley       maxRowLen = PetscMax(maxRowLen, len);
474*cb1e1211SMatthew G Knepley     }
475*cb1e1211SMatthew G Knepley     ierr = PetscMalloc(maxRowLen * sizeof(PetscScalar), &values);CHKERRQ(ierr);
476*cb1e1211SMatthew G Knepley     ierr = PetscMemzero(values, maxRowLen * sizeof(PetscScalar));CHKERRQ(ierr);
477*cb1e1211SMatthew G Knepley     for (r = rStart; r < rEnd; ++r) {
478*cb1e1211SMatthew G Knepley       PetscInt numCols, cStart;
479*cb1e1211SMatthew G Knepley 
480*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr);
481*cb1e1211SMatthew G Knepley       ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr);
482*cb1e1211SMatthew G Knepley       ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr);
483*cb1e1211SMatthew G Knepley     }
484*cb1e1211SMatthew G Knepley     ierr = PetscFree(values);CHKERRQ(ierr);
485*cb1e1211SMatthew G Knepley     ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
486*cb1e1211SMatthew G Knepley     ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
487*cb1e1211SMatthew G Knepley   }
488*cb1e1211SMatthew G Knepley   ierr = PetscSectionDestroy(&sectionAdj);CHKERRQ(ierr);
489*cb1e1211SMatthew G Knepley   ierr = PetscFree(cols);CHKERRQ(ierr);
490*cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
491*cb1e1211SMatthew G Knepley }
492*cb1e1211SMatthew G Knepley 
493*cb1e1211SMatthew G Knepley #if 0
494*cb1e1211SMatthew G Knepley #undef __FUNCT__
495*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator_2"
496*cb1e1211SMatthew 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)
497*cb1e1211SMatthew G Knepley {
498*cb1e1211SMatthew G Knepley   PetscInt       *tmpClosure,*tmpAdj,*visits;
499*cb1e1211SMatthew G Knepley   PetscInt        c,cStart,cEnd,pStart,pEnd;
500*cb1e1211SMatthew G Knepley   PetscErrorCode  ierr;
501*cb1e1211SMatthew G Knepley 
502*cb1e1211SMatthew G Knepley   PetscFunctionBegin;
503*cb1e1211SMatthew G Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
504*cb1e1211SMatthew G Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
505*cb1e1211SMatthew G Knepley   ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
506*cb1e1211SMatthew G Knepley 
507*cb1e1211SMatthew G Knepley   maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth),PetscPowInt(mesh->maxSupportSize,depth));
508*cb1e1211SMatthew G Knepley 
509*cb1e1211SMatthew G Knepley   ierr    = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr);
510*cb1e1211SMatthew G Knepley   npoints = pEnd - pStart;
511*cb1e1211SMatthew G Knepley 
512*cb1e1211SMatthew G Knepley   ierr = PetscMalloc3(maxClosureSize,PetscInt,&tmpClosure,npoints,PetscInt,&lvisits,npoints,PetscInt,&visits);CHKERRQ(ierr);
513*cb1e1211SMatthew G Knepley   ierr = PetscMemzero(lvisits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr);
514*cb1e1211SMatthew G Knepley   ierr = PetscMemzero(visits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr);
515*cb1e1211SMatthew G Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
516*cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
517*cb1e1211SMatthew G Knepley     PetscInt *support = tmpClosure;
518*cb1e1211SMatthew G Knepley     ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support);CHKERRQ(ierr);
519*cb1e1211SMatthew G Knepley     for (p=0; p<supportSize; p++) lvisits[support[p]]++;
520*cb1e1211SMatthew G Knepley   }
521*cb1e1211SMatthew G Knepley   ierr = PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr);
522*cb1e1211SMatthew G Knepley   ierr = PetscSFReduceEnd  (sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr);
523*cb1e1211SMatthew G Knepley   ierr = PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr);
524*cb1e1211SMatthew G Knepley   ierr = PetscSFBcastEnd  (sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr);
525*cb1e1211SMatthew G Knepley 
526*cb1e1211SMatthew G Knepley   ierr = PetscSFGetRanks();CHKERRQ(ierr);
527*cb1e1211SMatthew G Knepley 
528*cb1e1211SMatthew G Knepley 
529*cb1e1211SMatthew G Knepley   ierr = PetscMalloc2(maxClosureSize*maxClosureSize,PetscInt,&cellmat,npoints,PetscInt,&owner);CHKERRQ(ierr);
530*cb1e1211SMatthew G Knepley   for (c=cStart; c<cEnd; c++) {
531*cb1e1211SMatthew G Knepley     ierr = PetscMemzero(cellmat,maxClosureSize*maxClosureSize*sizeof(PetscInt));CHKERRQ(ierr);
532*cb1e1211SMatthew G Knepley     /*
533*cb1e1211SMatthew G Knepley      Depth-first walk of transitive closure.
534*cb1e1211SMatthew 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.
535*cb1e1211SMatthew G Knepley      This contribution is added to dnz if owning ranks of p and q match, to onz otherwise.
536*cb1e1211SMatthew G Knepley      */
537*cb1e1211SMatthew G Knepley   }
538*cb1e1211SMatthew G Knepley 
539*cb1e1211SMatthew G Knepley   ierr = PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM);CHKERRQ(ierr);
540*cb1e1211SMatthew G Knepley   ierr = PetscSFReduceEnd  (sf,MPIU_INT,lonz,onz,MPI_SUM);CHKERRQ(ierr);
541*cb1e1211SMatthew G Knepley   PetscFunctionReturn(0);
542*cb1e1211SMatthew G Knepley }
543*cb1e1211SMatthew G Knepley #endif
544