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__ "DMPlexPreallocateOperator" 7cb1e1211SMatthew G Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix) 8cb1e1211SMatthew G Knepley { 9cb1e1211SMatthew G Knepley MPI_Comm comm; 1089545effSMatthew G. Knepley MatType mtype; 11cb1e1211SMatthew G Knepley PetscSF sf, sfDof, sfAdj; 12cb1e1211SMatthew G Knepley PetscSection leafSectionAdj, rootSectionAdj, sectionAdj; 1347a6104aSMatthew G. Knepley PetscInt nroots, nleaves, l, p; 14cb1e1211SMatthew G Knepley const PetscInt *leaves; 15cb1e1211SMatthew G Knepley const PetscSFNode *remotes; 16cb1e1211SMatthew G Knepley PetscInt dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols; 1770034214SMatthew G. Knepley PetscInt *tmpAdj = NULL, *adj, *rootAdj, *cols, *remoteOffsets; 1870034214SMatthew G. Knepley PetscInt adjSize; 19cb1e1211SMatthew G Knepley PetscLayout rLayout; 20cb1e1211SMatthew G Knepley PetscInt locRows, rStart, rEnd, r; 21cb1e1211SMatthew G Knepley PetscMPIInt size; 2270034214SMatthew G. Knepley PetscBool doCommLocal, doComm, debug = PETSC_FALSE, isSymBlock, isSymSeqBlock, isSymMPIBlock; 23cb1e1211SMatthew G Knepley PetscErrorCode ierr; 24cb1e1211SMatthew G Knepley 25cb1e1211SMatthew G Knepley PetscFunctionBegin; 2670034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2770034214SMatthew G. Knepley PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 2870034214SMatthew G. Knepley PetscValidHeaderSpecific(sectionGlobal, PETSC_SECTION_CLASSID, 4); 2970034214SMatthew G. Knepley PetscValidHeaderSpecific(A, MAT_CLASSID, 9); 308c598124SMatthew G. Knepley if (dnz) PetscValidPointer(dnz,5); 318c598124SMatthew G. Knepley if (onz) PetscValidPointer(onz,6); 328c598124SMatthew G. Knepley if (dnzu) PetscValidPointer(dnzu,7); 338c598124SMatthew G. Knepley if (onzu) PetscValidPointer(onzu,8); 34a72f3261SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr); 35cb1e1211SMatthew G Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 36cb1e1211SMatthew G Knepley ierr = PetscOptionsGetBool(NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr); 37cb1e1211SMatthew G Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 38*c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 39cb1e1211SMatthew G Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 4047a6104aSMatthew G. Knepley ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 4147a6104aSMatthew G. Knepley doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE; 4247a6104aSMatthew G. Knepley ierr = MPI_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr); 43cb1e1211SMatthew G Knepley /* Create dof SF based on point SF */ 44cb1e1211SMatthew G Knepley if (debug) { 45cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Input Section for Preallocation:\n");CHKERRQ(ierr); 46cb1e1211SMatthew G Knepley ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 47cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Input Global Section for Preallocation:\n");CHKERRQ(ierr); 48cb1e1211SMatthew G Knepley ierr = PetscSectionView(sectionGlobal, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 49cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Input SF for Preallocation:\n");CHKERRQ(ierr); 50cb1e1211SMatthew G Knepley ierr = PetscSFView(sf, NULL);CHKERRQ(ierr); 51cb1e1211SMatthew G Knepley } 52cb1e1211SMatthew G Knepley ierr = PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets);CHKERRQ(ierr); 53cb1e1211SMatthew G Knepley ierr = PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof);CHKERRQ(ierr); 54cb1e1211SMatthew G Knepley if (debug) { 55cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Dof SF for Preallocation:\n");CHKERRQ(ierr); 56cb1e1211SMatthew G Knepley ierr = PetscSFView(sfDof, NULL);CHKERRQ(ierr); 57cb1e1211SMatthew G Knepley } 58cb1e1211SMatthew G Knepley /* Create section for dof adjacency (dof ==> # adj dof) */ 59cb1e1211SMatthew G Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 60cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(section, &numDof);CHKERRQ(ierr); 61cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, &leafSectionAdj);CHKERRQ(ierr); 62cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(leafSectionAdj, 0, numDof);CHKERRQ(ierr); 63cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, &rootSectionAdj);CHKERRQ(ierr); 64cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(rootSectionAdj, 0, numDof);CHKERRQ(ierr); 65cb1e1211SMatthew G Knepley /* Fill in the ghost dofs on the interface */ 66cb1e1211SMatthew G Knepley ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes);CHKERRQ(ierr); 67cb1e1211SMatthew G Knepley 68cb1e1211SMatthew G Knepley /* 69964bf7afSMatthew G. Knepley section - maps points to (# dofs, local dofs) 70964bf7afSMatthew G. Knepley sectionGlobal - maps points to (# dofs, global dofs) 71964bf7afSMatthew G. Knepley leafSectionAdj - maps unowned local dofs to # adj dofs 72964bf7afSMatthew G. Knepley rootSectionAdj - maps owned local dofs to # adj dofs 73964bf7afSMatthew G. Knepley adj - adj global dofs indexed by leafSectionAdj 74964bf7afSMatthew G. Knepley rootAdj - adj global dofs indexed by rootSectionAdj 75964bf7afSMatthew G. Knepley sf - describes shared points across procs 76964bf7afSMatthew G. Knepley sfDof - describes shared dofs across procs 77964bf7afSMatthew G. Knepley sfAdj - describes shared adjacent dofs across procs 78cb1e1211SMatthew G Knepley ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point. 79cb1e1211SMatthew G Knepley 1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj 80cb1e1211SMatthew G Knepley Reduce those counts to rootSectionAdj (now redundantly counting some interface points) 81cb1e1211SMatthew G Knepley 2. Visit owned points on interface, count adjacencies placing in rootSectionAdj 82cb1e1211SMatthew G Knepley Create sfAdj connecting rootSectionAdj and leafSectionAdj 83cb1e1211SMatthew G Knepley 3. Visit unowned points on interface, write adjacencies to adj 84cb1e1211SMatthew G Knepley Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies) 85cb1e1211SMatthew G Knepley 4. Visit owned points on interface, write adjacencies to rootAdj 86cb1e1211SMatthew G Knepley Remove redundancy in rootAdj 87cb1e1211SMatthew G Knepley ** The last two traversals use transitive closure 88cb1e1211SMatthew G Knepley 5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj) 89cb1e1211SMatthew G Knepley Allocate memory addressed by sectionAdj (cols) 90cb1e1211SMatthew G Knepley 6. Visit all owned points in the subdomain, insert dof adjacencies into cols 91cb1e1211SMatthew G Knepley ** Knowing all the column adjacencies, check ownership and sum into dnz and onz 92cb1e1211SMatthew G Knepley */ 93cb1e1211SMatthew G Knepley 94cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 95cb1e1211SMatthew G Knepley PetscInt dof, off, d, q; 9670034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 97cb1e1211SMatthew G Knepley 98fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 99cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 100cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 10170034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 102cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 103cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 104cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 105cb1e1211SMatthew G Knepley 106cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 107cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 108cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 109cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 110cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 111cb1e1211SMatthew G Knepley } 112cb1e1211SMatthew G Knepley } 113cb1e1211SMatthew G Knepley } 114cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr); 115cb1e1211SMatthew G Knepley if (debug) { 116cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n");CHKERRQ(ierr); 117cb1e1211SMatthew G Knepley ierr = PetscSectionView(leafSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 118cb1e1211SMatthew G Knepley } 119cb1e1211SMatthew G Knepley /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */ 12047a6104aSMatthew G. Knepley if (doComm) { 121cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr); 122cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr); 123cb1e1211SMatthew G Knepley } 124cb1e1211SMatthew G Knepley if (debug) { 125cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n");CHKERRQ(ierr); 126cb1e1211SMatthew G Knepley ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 127cb1e1211SMatthew G Knepley } 128cb1e1211SMatthew G Knepley /* Add in local adjacency sizes for owned dofs on interface (roots) */ 129cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 13070034214SMatthew G. Knepley PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q; 131cb1e1211SMatthew G Knepley 132cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 133cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 134cb1e1211SMatthew G Knepley if (!dof) continue; 135cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 136cb1e1211SMatthew G Knepley if (adof <= 0) continue; 13770034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 138cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 139cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 140cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 141cb1e1211SMatthew G Knepley 142cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 143cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 144cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 145cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 146cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 147cb1e1211SMatthew G Knepley } 148cb1e1211SMatthew G Knepley } 149cb1e1211SMatthew G Knepley } 150cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr); 151cb1e1211SMatthew G Knepley if (debug) { 152cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n");CHKERRQ(ierr); 153cb1e1211SMatthew G Knepley ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 154cb1e1211SMatthew G Knepley } 155cb1e1211SMatthew G Knepley /* Create adj SF based on dof SF */ 156cb1e1211SMatthew G Knepley ierr = PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets);CHKERRQ(ierr); 157cb1e1211SMatthew G Knepley ierr = PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj);CHKERRQ(ierr); 158cb1e1211SMatthew G Knepley if (debug) { 159cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency SF for Preallocation:\n");CHKERRQ(ierr); 160cb1e1211SMatthew G Knepley ierr = PetscSFView(sfAdj, NULL);CHKERRQ(ierr); 161cb1e1211SMatthew G Knepley } 162cb1e1211SMatthew G Knepley ierr = PetscSFDestroy(&sfDof);CHKERRQ(ierr); 163cb1e1211SMatthew G Knepley /* Create leaf adjacency */ 164cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr); 165cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(leafSectionAdj, &adjSize);CHKERRQ(ierr); 1661795a4d1SJed Brown ierr = PetscCalloc1(adjSize, &adj);CHKERRQ(ierr); 167cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 168cb1e1211SMatthew G Knepley PetscInt dof, off, d, q; 16970034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 170cb1e1211SMatthew G Knepley 171fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 172cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 173cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 17470034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 175cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 176cb1e1211SMatthew G Knepley PetscInt aoff, i = 0; 177cb1e1211SMatthew G Knepley 178cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(leafSectionAdj, d, &aoff);CHKERRQ(ierr); 179cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 180cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 181cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 182cb1e1211SMatthew G Knepley 183cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 184cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 185cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 186cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 187cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd) { 188cb1e1211SMatthew G Knepley adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd; 189cb1e1211SMatthew G Knepley ++i; 190cb1e1211SMatthew G Knepley } 191cb1e1211SMatthew G Knepley } 192cb1e1211SMatthew G Knepley } 193cb1e1211SMatthew G Knepley } 194cb1e1211SMatthew G Knepley /* Debugging */ 195cb1e1211SMatthew G Knepley if (debug) { 196cb1e1211SMatthew G Knepley IS tmp; 197cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Leaf adjacency indices\n");CHKERRQ(ierr); 198cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 199cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 2000a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 201cb1e1211SMatthew G Knepley } 202cb1e1211SMatthew G Knepley /* Gather adjacenct indices to root */ 203cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(rootSectionAdj, &adjSize);CHKERRQ(ierr); 204785e854fSJed Brown ierr = PetscMalloc1(adjSize, &rootAdj);CHKERRQ(ierr); 205cb1e1211SMatthew G Knepley for (r = 0; r < adjSize; ++r) rootAdj[r] = -1; 20647a6104aSMatthew G. Knepley if (doComm) { 207cb1e1211SMatthew G Knepley ierr = PetscSFGatherBegin(sfAdj, MPIU_INT, adj, rootAdj);CHKERRQ(ierr); 208cb1e1211SMatthew G Knepley ierr = PetscSFGatherEnd(sfAdj, MPIU_INT, adj, rootAdj);CHKERRQ(ierr); 209cb1e1211SMatthew G Knepley } 210cb1e1211SMatthew G Knepley ierr = PetscSFDestroy(&sfAdj);CHKERRQ(ierr); 211cb1e1211SMatthew G Knepley ierr = PetscFree(adj);CHKERRQ(ierr); 212cb1e1211SMatthew G Knepley /* Debugging */ 213cb1e1211SMatthew G Knepley if (debug) { 214cb1e1211SMatthew G Knepley IS tmp; 215cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices after gather\n");CHKERRQ(ierr); 216cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 217cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 2180a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 219cb1e1211SMatthew G Knepley } 220cb1e1211SMatthew G Knepley /* Add in local adjacency indices for owned dofs on interface (roots) */ 221cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 22270034214SMatthew G. Knepley PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q; 223cb1e1211SMatthew G Knepley 224cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 225cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 226cb1e1211SMatthew G Knepley if (!dof) continue; 227cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 228cb1e1211SMatthew G Knepley if (adof <= 0) continue; 22970034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 230cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 231cb1e1211SMatthew G Knepley PetscInt adof, aoff, i; 232cb1e1211SMatthew G Knepley 233cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr); 234cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr); 235cb1e1211SMatthew G Knepley i = adof-1; 236cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 237cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 238cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 239cb1e1211SMatthew G Knepley 240cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 241cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 242cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 243cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 244cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd) { 245cb1e1211SMatthew G Knepley rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd; 246cb1e1211SMatthew G Knepley --i; 247cb1e1211SMatthew G Knepley } 248cb1e1211SMatthew G Knepley } 249cb1e1211SMatthew G Knepley } 250cb1e1211SMatthew G Knepley } 251cb1e1211SMatthew G Knepley /* Debugging */ 252cb1e1211SMatthew G Knepley if (debug) { 253cb1e1211SMatthew G Knepley IS tmp; 254cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices\n");CHKERRQ(ierr); 255cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 256cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 2570a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 258cb1e1211SMatthew G Knepley } 259cb1e1211SMatthew G Knepley /* Compress indices */ 260cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr); 261cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 262cb1e1211SMatthew G Knepley PetscInt dof, cdof, off, d; 263cb1e1211SMatthew G Knepley PetscInt adof, aoff; 264cb1e1211SMatthew G Knepley 265cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 266cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 267cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 268cb1e1211SMatthew G Knepley if (!dof) continue; 269cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 270cb1e1211SMatthew G Knepley if (adof <= 0) continue; 271cb1e1211SMatthew G Knepley for (d = off; d < off+dof-cdof; ++d) { 272cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr); 273cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr); 274cb1e1211SMatthew G Knepley ierr = PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]);CHKERRQ(ierr); 275cb1e1211SMatthew G Knepley ierr = PetscSectionSetDof(rootSectionAdj, d, adof);CHKERRQ(ierr); 276cb1e1211SMatthew G Knepley } 277cb1e1211SMatthew G Knepley } 278cb1e1211SMatthew G Knepley /* Debugging */ 279cb1e1211SMatthew G Knepley if (debug) { 280cb1e1211SMatthew G Knepley IS tmp; 281cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n");CHKERRQ(ierr); 282cb1e1211SMatthew G Knepley ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 283cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices after compression\n");CHKERRQ(ierr); 284cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 285cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 2860a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 287cb1e1211SMatthew G Knepley } 288cb1e1211SMatthew G Knepley /* Build adjacency section: Maps global indices to sets of adjacent global indices */ 289cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd);CHKERRQ(ierr); 290cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, §ionAdj);CHKERRQ(ierr); 291cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd);CHKERRQ(ierr); 292cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 29370034214SMatthew G. Knepley PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q; 294cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 295cb1e1211SMatthew G Knepley 296cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 297cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 298cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 299cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 300cb1e1211SMatthew G Knepley for (d = 0; d < dof-cdof; ++d) { 301cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 302cb1e1211SMatthew G Knepley 303cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr); 304cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr); 305cb1e1211SMatthew G Knepley if (ldof > 0) { 306cb1e1211SMatthew G Knepley /* We do not own this point */ 307cb1e1211SMatthew G Knepley } else if (rdof > 0) { 308cb1e1211SMatthew G Knepley ierr = PetscSectionSetDof(sectionAdj, goff+d, rdof);CHKERRQ(ierr); 309cb1e1211SMatthew G Knepley } else { 310cb1e1211SMatthew G Knepley found = PETSC_FALSE; 311cb1e1211SMatthew G Knepley } 312cb1e1211SMatthew G Knepley } 313cb1e1211SMatthew G Knepley if (found) continue; 314cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 315cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 31670034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 317cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 318cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 319cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, noff; 320cb1e1211SMatthew G Knepley 321cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 322cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 323cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 324cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, padj, &noff);CHKERRQ(ierr); 325cb1e1211SMatthew G Knepley for (d = goff; d < goff+dof-cdof; ++d) { 326cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(sectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 327cb1e1211SMatthew G Knepley } 328cb1e1211SMatthew G Knepley } 329cb1e1211SMatthew G Knepley } 330cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(sectionAdj);CHKERRQ(ierr); 331cb1e1211SMatthew G Knepley if (debug) { 332cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency Section for Preallocation:\n");CHKERRQ(ierr); 333cb1e1211SMatthew G Knepley ierr = PetscSectionView(sectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 334cb1e1211SMatthew G Knepley } 335cb1e1211SMatthew G Knepley /* Get adjacent indices */ 336cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(sectionAdj, &numCols);CHKERRQ(ierr); 337785e854fSJed Brown ierr = PetscMalloc1(numCols, &cols);CHKERRQ(ierr); 338cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 33970034214SMatthew G. Knepley PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q; 340cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 341cb1e1211SMatthew G Knepley 342cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 343cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 344cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 345cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 346cb1e1211SMatthew G Knepley for (d = 0; d < dof-cdof; ++d) { 347cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 348cb1e1211SMatthew G Knepley 349cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr); 350cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr); 351cb1e1211SMatthew G Knepley if (ldof > 0) { 352cb1e1211SMatthew G Knepley /* We do not own this point */ 353cb1e1211SMatthew G Knepley } else if (rdof > 0) { 354cb1e1211SMatthew G Knepley PetscInt aoff, roff; 355cb1e1211SMatthew G Knepley 356cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, goff+d, &aoff);CHKERRQ(ierr); 357cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, off+d, &roff);CHKERRQ(ierr); 358cb1e1211SMatthew G Knepley ierr = PetscMemcpy(&cols[aoff], &rootAdj[roff], rdof * sizeof(PetscInt));CHKERRQ(ierr); 359cb1e1211SMatthew G Knepley } else { 360cb1e1211SMatthew G Knepley found = PETSC_FALSE; 361cb1e1211SMatthew G Knepley } 362cb1e1211SMatthew G Knepley } 363cb1e1211SMatthew G Knepley if (found) continue; 36470034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 365cb1e1211SMatthew G Knepley for (d = goff; d < goff+dof-cdof; ++d) { 366cb1e1211SMatthew G Knepley PetscInt adof, aoff, i = 0; 367cb1e1211SMatthew G Knepley 368cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, d, &adof);CHKERRQ(ierr); 369cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, d, &aoff);CHKERRQ(ierr); 370cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 371cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 372cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 373cb1e1211SMatthew G Knepley const PetscInt *ncind; 374cb1e1211SMatthew G Knepley 375cb1e1211SMatthew G Knepley /* Adjacent points may not be in the section chart */ 376cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 377cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 378cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 379cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintIndices(section, padj, &ncind);CHKERRQ(ierr); 380cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 381cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd, ++i) { 382cb1e1211SMatthew G Knepley cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd; 383cb1e1211SMatthew G Knepley } 384cb1e1211SMatthew G Knepley } 385cb1e1211SMatthew 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); 386cb1e1211SMatthew G Knepley } 387cb1e1211SMatthew G Knepley } 388cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(&leafSectionAdj);CHKERRQ(ierr); 389cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(&rootSectionAdj);CHKERRQ(ierr); 390cb1e1211SMatthew G Knepley ierr = PetscFree(rootAdj);CHKERRQ(ierr); 39170034214SMatthew G. Knepley ierr = PetscFree(tmpAdj);CHKERRQ(ierr); 392cb1e1211SMatthew G Knepley /* Debugging */ 393cb1e1211SMatthew G Knepley if (debug) { 394cb1e1211SMatthew G Knepley IS tmp; 395cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Column indices\n");CHKERRQ(ierr); 396cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 397cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 3980a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 399cb1e1211SMatthew G Knepley } 400cb1e1211SMatthew G Knepley /* Create allocation vectors from adjacency graph */ 401cb1e1211SMatthew G Knepley ierr = MatGetLocalSize(A, &locRows, NULL);CHKERRQ(ierr); 402cb1e1211SMatthew G Knepley ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)A), &rLayout);CHKERRQ(ierr); 403cb1e1211SMatthew G Knepley ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 404cb1e1211SMatthew G Knepley ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 405cb1e1211SMatthew G Knepley ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 406cb1e1211SMatthew G Knepley ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 407cb1e1211SMatthew G Knepley ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 408cb1e1211SMatthew G Knepley /* Only loop over blocks of rows */ 409cb1e1211SMatthew 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); 410cb1e1211SMatthew G Knepley for (r = rStart/bs; r < rEnd/bs; ++r) { 411cb1e1211SMatthew G Knepley const PetscInt row = r*bs; 412cb1e1211SMatthew G Knepley PetscInt numCols, cStart, c; 413cb1e1211SMatthew G Knepley 414cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, row, &numCols);CHKERRQ(ierr); 415cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, row, &cStart);CHKERRQ(ierr); 416cb1e1211SMatthew G Knepley for (c = cStart; c < cStart+numCols; ++c) { 417cb1e1211SMatthew G Knepley if ((cols[c] >= rStart*bs) && (cols[c] < rEnd*bs)) { 418cb1e1211SMatthew G Knepley ++dnz[r-rStart]; 419cb1e1211SMatthew G Knepley if (cols[c] >= row) ++dnzu[r-rStart]; 420cb1e1211SMatthew G Knepley } else { 421cb1e1211SMatthew G Knepley ++onz[r-rStart]; 422cb1e1211SMatthew G Knepley if (cols[c] >= row) ++onzu[r-rStart]; 423cb1e1211SMatthew G Knepley } 424cb1e1211SMatthew G Knepley } 425cb1e1211SMatthew G Knepley } 426cb1e1211SMatthew G Knepley if (bs > 1) { 427cb1e1211SMatthew G Knepley for (r = 0; r < locRows/bs; ++r) { 428cb1e1211SMatthew G Knepley dnz[r] /= bs; 429cb1e1211SMatthew G Knepley onz[r] /= bs; 430cb1e1211SMatthew G Knepley dnzu[r] /= bs; 431cb1e1211SMatthew G Knepley onzu[r] /= bs; 432cb1e1211SMatthew G Knepley } 433cb1e1211SMatthew G Knepley } 434cb1e1211SMatthew G Knepley /* Set matrix pattern */ 435cb1e1211SMatthew G Knepley ierr = MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu);CHKERRQ(ierr); 436cb1e1211SMatthew G Knepley ierr = MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 43789545effSMatthew G. Knepley /* Check for symmetric storage */ 43889545effSMatthew G. Knepley ierr = MatGetType(A, &mtype);CHKERRQ(ierr); 43989545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr); 44089545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr); 44189545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr); 44289545effSMatthew G. Knepley if (isSymBlock || isSymSeqBlock || isSymMPIBlock) {ierr = MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);} 443cb1e1211SMatthew G Knepley /* Fill matrix with zeros */ 444cb1e1211SMatthew G Knepley if (fillMatrix) { 445cb1e1211SMatthew G Knepley PetscScalar *values; 446cb1e1211SMatthew G Knepley PetscInt maxRowLen = 0; 447cb1e1211SMatthew G Knepley 448cb1e1211SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 449cb1e1211SMatthew G Knepley PetscInt len; 450cb1e1211SMatthew G Knepley 451cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, r, &len);CHKERRQ(ierr); 452cb1e1211SMatthew G Knepley maxRowLen = PetscMax(maxRowLen, len); 453cb1e1211SMatthew G Knepley } 4541795a4d1SJed Brown ierr = PetscCalloc1(maxRowLen, &values);CHKERRQ(ierr); 455cb1e1211SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 456cb1e1211SMatthew G Knepley PetscInt numCols, cStart; 457cb1e1211SMatthew G Knepley 458cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr); 459cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr); 460cb1e1211SMatthew G Knepley ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr); 461cb1e1211SMatthew G Knepley } 462cb1e1211SMatthew G Knepley ierr = PetscFree(values);CHKERRQ(ierr); 463cb1e1211SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 464cb1e1211SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 465cb1e1211SMatthew G Knepley } 466cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(§ionAdj);CHKERRQ(ierr); 467cb1e1211SMatthew G Knepley ierr = PetscFree(cols);CHKERRQ(ierr); 468a72f3261SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr); 469cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 470cb1e1211SMatthew G Knepley } 471cb1e1211SMatthew G Knepley 472cb1e1211SMatthew G Knepley #if 0 473cb1e1211SMatthew G Knepley #undef __FUNCT__ 474cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator_2" 475cb1e1211SMatthew 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) 476cb1e1211SMatthew G Knepley { 477cb1e1211SMatthew G Knepley PetscInt *tmpClosure,*tmpAdj,*visits; 478cb1e1211SMatthew G Knepley PetscInt c,cStart,cEnd,pStart,pEnd; 479cb1e1211SMatthew G Knepley PetscErrorCode ierr; 480cb1e1211SMatthew G Knepley 481cb1e1211SMatthew G Knepley PetscFunctionBegin; 482*c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 483cb1e1211SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 484cb1e1211SMatthew G Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 485cb1e1211SMatthew G Knepley 486e7b80042SMatthew G. Knepley maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1)); 487cb1e1211SMatthew G Knepley 488cb1e1211SMatthew G Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 489cb1e1211SMatthew G Knepley npoints = pEnd - pStart; 490cb1e1211SMatthew G Knepley 491dcca6d9dSJed Brown ierr = PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits);CHKERRQ(ierr); 492cb1e1211SMatthew G Knepley ierr = PetscMemzero(lvisits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr); 493cb1e1211SMatthew G Knepley ierr = PetscMemzero(visits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr); 494cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 495cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 496cb1e1211SMatthew G Knepley PetscInt *support = tmpClosure; 497cb1e1211SMatthew G Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support);CHKERRQ(ierr); 498cb1e1211SMatthew G Knepley for (p=0; p<supportSize; p++) lvisits[support[p]]++; 499cb1e1211SMatthew G Knepley } 500cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr); 501cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd (sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr); 502cb1e1211SMatthew G Knepley ierr = PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr); 503cb1e1211SMatthew G Knepley ierr = PetscSFBcastEnd (sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr); 504cb1e1211SMatthew G Knepley 505cb1e1211SMatthew G Knepley ierr = PetscSFGetRanks();CHKERRQ(ierr); 506cb1e1211SMatthew G Knepley 507cb1e1211SMatthew G Knepley 508dcca6d9dSJed Brown ierr = PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner);CHKERRQ(ierr); 509cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 510cb1e1211SMatthew G Knepley ierr = PetscMemzero(cellmat,maxClosureSize*maxClosureSize*sizeof(PetscInt));CHKERRQ(ierr); 511cb1e1211SMatthew G Knepley /* 512cb1e1211SMatthew G Knepley Depth-first walk of transitive closure. 513cb1e1211SMatthew 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. 514cb1e1211SMatthew G Knepley This contribution is added to dnz if owning ranks of p and q match, to onz otherwise. 515cb1e1211SMatthew G Knepley */ 516cb1e1211SMatthew G Knepley } 517cb1e1211SMatthew G Knepley 518cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM);CHKERRQ(ierr); 519cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd (sf,MPIU_INT,lonz,onz,MPI_SUM);CHKERRQ(ierr); 520cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 521cb1e1211SMatthew G Knepley } 522cb1e1211SMatthew G Knepley #endif 523