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