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