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