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__ 6*76185916SToby Isaac #define __FUNCT__ "DMPlexComputeAnchorAdjacencies" 7*76185916SToby Isaac /* get adjacencies due to point-to-point constraints that can't be found with DMPlexGetAdjacency() */ 8*76185916SToby Isaac static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscSection section, PetscSection sectionGlobal, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[]) 9*76185916SToby Isaac { 10*76185916SToby Isaac PetscInt pStart, pEnd; 11*76185916SToby Isaac PetscSection adjSec, aSec; 12*76185916SToby Isaac IS aIS; 13*76185916SToby Isaac PetscErrorCode ierr; 14*76185916SToby Isaac 15*76185916SToby Isaac PetscFunctionBegin; 16*76185916SToby Isaac 17*76185916SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&adjSec);CHKERRQ(ierr); 18*76185916SToby Isaac ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 19*76185916SToby Isaac ierr = PetscSectionSetChart(adjSec,pStart,pEnd);CHKERRQ(ierr); 20*76185916SToby Isaac 21*76185916SToby Isaac ierr = DMPlexGetConstraints(dm,&aSec,&aIS);CHKERRQ(ierr); 22*76185916SToby Isaac if (aSec) { 23*76185916SToby Isaac const PetscInt *anchors; 24*76185916SToby Isaac PetscInt p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize; 25*76185916SToby Isaac PetscInt *tmpAdjP = NULL, *tmpAdjQ = NULL; 26*76185916SToby Isaac PetscSection inverseSec; 27*76185916SToby Isaac 28*76185916SToby Isaac /* invert the constraint-to-anchor map */ 29*76185916SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)aSec),&inverseSec);CHKERRQ(ierr); 30*76185916SToby Isaac ierr = PetscSectionSetChart(inverseSec,pStart,pEnd);CHKERRQ(ierr); 31*76185916SToby Isaac ierr = ISGetLocalSize(aIS, &aSize);CHKERRQ(ierr); 32*76185916SToby Isaac ierr = ISGetIndices(aIS, &anchors);CHKERRQ(ierr); 33*76185916SToby Isaac 34*76185916SToby Isaac for (p = 0; p < aSize; p++) { 35*76185916SToby Isaac PetscInt a = anchors[p]; 36*76185916SToby Isaac 37*76185916SToby Isaac ierr = PetscSectionAddDof(inverseSec,a,1);CHKERRQ(ierr); 38*76185916SToby Isaac } 39*76185916SToby Isaac ierr = PetscSectionSetUp(inverseSec);CHKERRQ(ierr); 40*76185916SToby Isaac ierr = PetscSectionGetStorageSize(inverseSec,&iSize);CHKERRQ(ierr); 41*76185916SToby Isaac ierr = PetscMalloc1(iSize,&inverse);CHKERRQ(ierr); 42*76185916SToby Isaac ierr = PetscCalloc1(pEnd-pStart,&offsets);CHKERRQ(ierr); 43*76185916SToby Isaac ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr); 44*76185916SToby Isaac for (p = aStart; p < aEnd; p++) { 45*76185916SToby Isaac PetscInt dof, off; 46*76185916SToby Isaac 47*76185916SToby Isaac ierr = PetscSectionGetDof(aSec, p, &dof);CHKERRQ(ierr); 48*76185916SToby Isaac ierr = PetscSectionGetOffset(aSec, p, &off);CHKERRQ(ierr); 49*76185916SToby Isaac 50*76185916SToby Isaac for (q = 0; q < dof; q++) { 51*76185916SToby Isaac PetscInt iOff; 52*76185916SToby Isaac 53*76185916SToby Isaac a = anchors[off + q]; 54*76185916SToby Isaac ierr = PetscSectionGetOffset(inverseSec, a, &iOff);CHKERRQ(ierr); 55*76185916SToby Isaac inverse[iOff + offsets[a-pStart]++] = p; 56*76185916SToby Isaac } 57*76185916SToby Isaac } 58*76185916SToby Isaac ierr = ISRestoreIndices(aIS, &anchors);CHKERRQ(ierr); 59*76185916SToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 60*76185916SToby Isaac 61*76185916SToby Isaac /* construct anchorAdj and adjSec 62*76185916SToby Isaac * 63*76185916SToby Isaac * loop over anchors: 64*76185916SToby Isaac * construct anchor adjacency 65*76185916SToby Isaac * loop over constrained: 66*76185916SToby Isaac * construct constrained adjacency 67*76185916SToby Isaac * if not in anchor adjacency, add to dofs 68*76185916SToby Isaac * setup adjSec, allocate anchorAdj 69*76185916SToby Isaac * loop over anchors: 70*76185916SToby Isaac * construct anchor adjacency 71*76185916SToby Isaac * loop over constrained: 72*76185916SToby Isaac * construct constrained adjacency 73*76185916SToby Isaac * if not in anchor adjacency 74*76185916SToby Isaac * if not already in list, put in list 75*76185916SToby Isaac * sort, unique, reduce dof count 76*76185916SToby Isaac * optional: compactify 77*76185916SToby Isaac */ 78*76185916SToby Isaac for (p = pStart; p < pEnd; p++) { 79*76185916SToby Isaac PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE; 80*76185916SToby Isaac 81*76185916SToby Isaac ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr); 82*76185916SToby Isaac if (!iDof) continue; 83*76185916SToby Isaac ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr); 84*76185916SToby Isaac ierr = DMPlexGetAdjacency(dm,p,&numAdjP,&tmpAdjP);CHKERRQ(ierr); 85*76185916SToby Isaac for (i = 0; i < iDof; i++) { 86*76185916SToby Isaac PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE; 87*76185916SToby Isaac 88*76185916SToby Isaac q = inverse[iOff + i]; 89*76185916SToby Isaac ierr = DMPlexGetAdjacency(dm,q,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr); 90*76185916SToby Isaac for (r = 0; r < numAdjQ; r++) { 91*76185916SToby Isaac qAdj = tmpAdjQ[r]; 92*76185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue; 93*76185916SToby Isaac for (s = 0; s < numAdjP; s++) { 94*76185916SToby Isaac if (qAdj == tmpAdjP[s]) break; 95*76185916SToby Isaac } 96*76185916SToby Isaac if (s < numAdjP) continue; 97*76185916SToby Isaac ierr = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr); 98*76185916SToby Isaac ierr = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr); 99*76185916SToby Isaac iNew += qAdjDof - qAdjCDof; 100*76185916SToby Isaac } 101*76185916SToby Isaac ierr = PetscSectionAddDof(adjSec,p,iNew);CHKERRQ(ierr); 102*76185916SToby Isaac } 103*76185916SToby Isaac } 104*76185916SToby Isaac 105*76185916SToby Isaac ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr); 106*76185916SToby Isaac ierr = PetscSectionGetStorageSize(adjSec,&adjSize);CHKERRQ(ierr); 107*76185916SToby Isaac ierr = PetscMalloc1(adjSize,&adj);CHKERRQ(ierr); 108*76185916SToby Isaac 109*76185916SToby Isaac for (p = pStart; p < pEnd; p++) { 110*76185916SToby Isaac PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE; 111*76185916SToby Isaac 112*76185916SToby Isaac ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr); 113*76185916SToby Isaac if (!iDof) continue; 114*76185916SToby Isaac ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr); 115*76185916SToby Isaac ierr = DMPlexGetAdjacency(dm,p,&numAdjP,&tmpAdjP);CHKERRQ(ierr); 116*76185916SToby Isaac ierr = PetscSectionGetDof(adjSec,p,&aDof);CHKERRQ(ierr); 117*76185916SToby Isaac ierr = PetscSectionGetOffset(adjSec,p,&aOff);CHKERRQ(ierr); 118*76185916SToby Isaac aOffOrig = aOff; 119*76185916SToby Isaac for (i = 0; i < iDof; i++) { 120*76185916SToby Isaac PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE; 121*76185916SToby Isaac 122*76185916SToby Isaac q = inverse[iOff + i]; 123*76185916SToby Isaac ierr = DMPlexGetAdjacency(dm,q,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr); 124*76185916SToby Isaac for (r = 0; r < numAdjQ; r++) { 125*76185916SToby Isaac qAdj = tmpAdjQ[r]; 126*76185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue; 127*76185916SToby Isaac for (s = 0; s < numAdjP; s++) { 128*76185916SToby Isaac if (qAdj == tmpAdjP[s]) break; 129*76185916SToby Isaac } 130*76185916SToby Isaac if (s < numAdjP) continue; 131*76185916SToby Isaac ierr = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr); 132*76185916SToby Isaac ierr = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr); 133*76185916SToby Isaac ierr = PetscSectionGetOffset(sectionGlobal,qAdj,&qAdjOff);CHKERRQ(ierr); 134*76185916SToby Isaac for (nd = 0; nd < qAdjDof-qAdjCDof; ++nd) { 135*76185916SToby Isaac adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff+1) : qAdjOff) + nd; 136*76185916SToby Isaac } 137*76185916SToby Isaac } 138*76185916SToby Isaac } 139*76185916SToby Isaac ierr = PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]);CHKERRQ(ierr); 140*76185916SToby Isaac ierr = PetscSectionSetDof(adjSec,p,aDof);CHKERRQ(ierr); 141*76185916SToby Isaac } 142*76185916SToby Isaac *anchorAdj = adj; 143*76185916SToby Isaac 144*76185916SToby Isaac /* clean up */ 145*76185916SToby Isaac ierr = PetscSectionDestroy(&inverseSec);CHKERRQ(ierr); 146*76185916SToby Isaac ierr = PetscFree(inverse);CHKERRQ(ierr); 147*76185916SToby Isaac ierr = PetscFree(tmpAdjP);CHKERRQ(ierr); 148*76185916SToby Isaac ierr = PetscFree(tmpAdjQ);CHKERRQ(ierr); 149*76185916SToby Isaac } 150*76185916SToby Isaac else { 151*76185916SToby Isaac *anchorAdj = NULL; 152*76185916SToby Isaac ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr); 153*76185916SToby Isaac } 154*76185916SToby Isaac *anchorSectionAdj = adjSec; 155*76185916SToby Isaac PetscFunctionReturn(0); 156*76185916SToby Isaac } 157*76185916SToby Isaac 158*76185916SToby Isaac #undef __FUNCT__ 159cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator" 160cb1e1211SMatthew G Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscSection section, PetscSection sectionGlobal, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix) 161cb1e1211SMatthew G Knepley { 162cb1e1211SMatthew G Knepley MPI_Comm comm; 16389545effSMatthew G. Knepley MatType mtype; 164cb1e1211SMatthew G Knepley PetscSF sf, sfDof, sfAdj; 165*76185916SToby Isaac PetscSection leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj; 16647a6104aSMatthew G. Knepley PetscInt nroots, nleaves, l, p; 167cb1e1211SMatthew G Knepley const PetscInt *leaves; 168cb1e1211SMatthew G Knepley const PetscSFNode *remotes; 169cb1e1211SMatthew G Knepley PetscInt dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols; 170*76185916SToby Isaac PetscInt *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj, *cols, *remoteOffsets; 17170034214SMatthew G. Knepley PetscInt adjSize; 172cb1e1211SMatthew G Knepley PetscLayout rLayout; 173cb1e1211SMatthew G Knepley PetscInt locRows, rStart, rEnd, r; 174cb1e1211SMatthew G Knepley PetscMPIInt size; 17570034214SMatthew G. Knepley PetscBool doCommLocal, doComm, debug = PETSC_FALSE, isSymBlock, isSymSeqBlock, isSymMPIBlock; 1768b0b4c70SToby Isaac PetscBool useConstraints; 177cb1e1211SMatthew G Knepley PetscErrorCode ierr; 178cb1e1211SMatthew G Knepley 179cb1e1211SMatthew G Knepley PetscFunctionBegin; 18070034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 18170034214SMatthew G. Knepley PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 3); 18270034214SMatthew G. Knepley PetscValidHeaderSpecific(sectionGlobal, PETSC_SECTION_CLASSID, 4); 18370034214SMatthew G. Knepley PetscValidHeaderSpecific(A, MAT_CLASSID, 9); 1848c598124SMatthew G. Knepley if (dnz) PetscValidPointer(dnz,5); 1858c598124SMatthew G. Knepley if (onz) PetscValidPointer(onz,6); 1868c598124SMatthew G. Knepley if (dnzu) PetscValidPointer(dnzu,7); 1878c598124SMatthew G. Knepley if (onzu) PetscValidPointer(onzu,8); 188a72f3261SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr); 189cb1e1211SMatthew G Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 190cb1e1211SMatthew G Knepley ierr = PetscOptionsGetBool(NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr); 191cb1e1211SMatthew G Knepley ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 192c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 193cb1e1211SMatthew G Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 19447a6104aSMatthew G. Knepley ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 19547a6104aSMatthew G. Knepley doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE; 19647a6104aSMatthew G. Knepley ierr = MPI_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRQ(ierr); 197cb1e1211SMatthew G Knepley /* Create dof SF based on point SF */ 198cb1e1211SMatthew G Knepley if (debug) { 199cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Input Section for Preallocation:\n");CHKERRQ(ierr); 200cb1e1211SMatthew G Knepley ierr = PetscSectionView(section, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 201cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Input Global Section for Preallocation:\n");CHKERRQ(ierr); 202cb1e1211SMatthew G Knepley ierr = PetscSectionView(sectionGlobal, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 203cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Input SF for Preallocation:\n");CHKERRQ(ierr); 204cb1e1211SMatthew G Knepley ierr = PetscSFView(sf, NULL);CHKERRQ(ierr); 205cb1e1211SMatthew G Knepley } 206cb1e1211SMatthew G Knepley ierr = PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets);CHKERRQ(ierr); 207cb1e1211SMatthew G Knepley ierr = PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof);CHKERRQ(ierr); 208cb1e1211SMatthew G Knepley if (debug) { 209cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Dof SF for Preallocation:\n");CHKERRQ(ierr); 210cb1e1211SMatthew G Knepley ierr = PetscSFView(sfDof, NULL);CHKERRQ(ierr); 211cb1e1211SMatthew G Knepley } 212cb1e1211SMatthew G Knepley /* Create section for dof adjacency (dof ==> # adj dof) */ 213cb1e1211SMatthew G Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 214cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(section, &numDof);CHKERRQ(ierr); 215cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, &leafSectionAdj);CHKERRQ(ierr); 216cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(leafSectionAdj, 0, numDof);CHKERRQ(ierr); 217cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, &rootSectionAdj);CHKERRQ(ierr); 218cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(rootSectionAdj, 0, numDof);CHKERRQ(ierr); 219cb1e1211SMatthew G Knepley /* Fill in the ghost dofs on the interface */ 220cb1e1211SMatthew G Knepley ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes);CHKERRQ(ierr); 2218b0b4c70SToby Isaac /* use constraints in finding adjacency in this routine */ 2228b0b4c70SToby Isaac ierr = DMPlexGetAdjacencyUseConstraints(dm,&useConstraints);CHKERRQ(ierr); 2238b0b4c70SToby Isaac ierr = DMPlexSetAdjacencyUseConstraints(dm,PETSC_TRUE);CHKERRQ(ierr); 224cb1e1211SMatthew G Knepley 225cb1e1211SMatthew G Knepley /* 226964bf7afSMatthew G. Knepley section - maps points to (# dofs, local dofs) 227964bf7afSMatthew G. Knepley sectionGlobal - maps points to (# dofs, global dofs) 228964bf7afSMatthew G. Knepley leafSectionAdj - maps unowned local dofs to # adj dofs 229964bf7afSMatthew G. Knepley rootSectionAdj - maps owned local dofs to # adj dofs 230964bf7afSMatthew G. Knepley adj - adj global dofs indexed by leafSectionAdj 231964bf7afSMatthew G. Knepley rootAdj - adj global dofs indexed by rootSectionAdj 232964bf7afSMatthew G. Knepley sf - describes shared points across procs 233964bf7afSMatthew G. Knepley sfDof - describes shared dofs across procs 234964bf7afSMatthew G. Knepley sfAdj - describes shared adjacent dofs across procs 235cb1e1211SMatthew G Knepley ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point. 236*76185916SToby Isaac (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj 237*76185916SToby Isaac (This is done in DMPlexComputeAnchorAdjacencies()) 238cb1e1211SMatthew G Knepley 1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj 239cb1e1211SMatthew G Knepley Reduce those counts to rootSectionAdj (now redundantly counting some interface points) 240cb1e1211SMatthew G Knepley 2. Visit owned points on interface, count adjacencies placing in rootSectionAdj 241cb1e1211SMatthew G Knepley Create sfAdj connecting rootSectionAdj and leafSectionAdj 242cb1e1211SMatthew G Knepley 3. Visit unowned points on interface, write adjacencies to adj 243cb1e1211SMatthew G Knepley Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies) 244cb1e1211SMatthew G Knepley 4. Visit owned points on interface, write adjacencies to rootAdj 245cb1e1211SMatthew G Knepley Remove redundancy in rootAdj 246cb1e1211SMatthew G Knepley ** The last two traversals use transitive closure 247cb1e1211SMatthew G Knepley 5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj) 248cb1e1211SMatthew G Knepley Allocate memory addressed by sectionAdj (cols) 249cb1e1211SMatthew G Knepley 6. Visit all owned points in the subdomain, insert dof adjacencies into cols 250cb1e1211SMatthew G Knepley ** Knowing all the column adjacencies, check ownership and sum into dnz and onz 251cb1e1211SMatthew G Knepley */ 252cb1e1211SMatthew G Knepley 253*76185916SToby Isaac ierr = DMPlexComputeAnchorAdjacencies(dm,section,sectionGlobal,&anchorSectionAdj,&anchorAdj);CHKERRQ(ierr); 254*76185916SToby Isaac 255cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 256*76185916SToby Isaac PetscInt dof, off, d, q, anDof; 25770034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 258cb1e1211SMatthew G Knepley 259fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 260cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 261cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 26270034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 263cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 264cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 265cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 266cb1e1211SMatthew G Knepley 267cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 268cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 269cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 270cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 271cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 272cb1e1211SMatthew G Knepley } 273cb1e1211SMatthew G Knepley } 274*76185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 275*76185916SToby Isaac if (anDof) { 276*76185916SToby Isaac for (d = off; d < off+dof; ++d) { 277*76185916SToby Isaac ierr = PetscSectionAddDof(leafSectionAdj, d, anDof);CHKERRQ(ierr); 278*76185916SToby Isaac } 279*76185916SToby Isaac } 280cb1e1211SMatthew G Knepley } 281cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr); 282cb1e1211SMatthew G Knepley if (debug) { 283cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n");CHKERRQ(ierr); 284cb1e1211SMatthew G Knepley ierr = PetscSectionView(leafSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 285cb1e1211SMatthew G Knepley } 286cb1e1211SMatthew G Knepley /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */ 28747a6104aSMatthew G. Knepley if (doComm) { 288cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr); 289cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr); 290cb1e1211SMatthew G Knepley } 291cb1e1211SMatthew G Knepley if (debug) { 292cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n");CHKERRQ(ierr); 293cb1e1211SMatthew G Knepley ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 294cb1e1211SMatthew G Knepley } 295cb1e1211SMatthew G Knepley /* Add in local adjacency sizes for owned dofs on interface (roots) */ 296cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 297*76185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof; 298cb1e1211SMatthew G Knepley 299cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 300cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 301cb1e1211SMatthew G Knepley if (!dof) continue; 302cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 303cb1e1211SMatthew G Knepley if (adof <= 0) continue; 30470034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 305cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 306cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 307cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 308cb1e1211SMatthew G Knepley 309cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 310cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 311cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 312cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 313cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 314cb1e1211SMatthew G Knepley } 315cb1e1211SMatthew G Knepley } 316*76185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 317*76185916SToby Isaac if (anDof) { 318*76185916SToby Isaac for (d = off; d < off+dof; ++d) { 319*76185916SToby Isaac ierr = PetscSectionAddDof(rootSectionAdj, d, anDof);CHKERRQ(ierr); 320*76185916SToby Isaac } 321*76185916SToby Isaac } 322cb1e1211SMatthew G Knepley } 323cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr); 324cb1e1211SMatthew G Knepley if (debug) { 325cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n");CHKERRQ(ierr); 326cb1e1211SMatthew G Knepley ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 327cb1e1211SMatthew G Knepley } 328cb1e1211SMatthew G Knepley /* Create adj SF based on dof SF */ 329cb1e1211SMatthew G Knepley ierr = PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets);CHKERRQ(ierr); 330cb1e1211SMatthew G Knepley ierr = PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj);CHKERRQ(ierr); 331cb1e1211SMatthew G Knepley if (debug) { 332cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency SF for Preallocation:\n");CHKERRQ(ierr); 333cb1e1211SMatthew G Knepley ierr = PetscSFView(sfAdj, NULL);CHKERRQ(ierr); 334cb1e1211SMatthew G Knepley } 335cb1e1211SMatthew G Knepley ierr = PetscSFDestroy(&sfDof);CHKERRQ(ierr); 336cb1e1211SMatthew G Knepley /* Create leaf adjacency */ 337cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr); 338cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(leafSectionAdj, &adjSize);CHKERRQ(ierr); 3391795a4d1SJed Brown ierr = PetscCalloc1(adjSize, &adj);CHKERRQ(ierr); 340cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 341*76185916SToby Isaac PetscInt dof, off, d, q, anDof, anOff; 34270034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 343cb1e1211SMatthew G Knepley 344fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 345cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 346cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 34770034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 348*76185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 349*76185916SToby Isaac ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr); 350cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 351cb1e1211SMatthew G Knepley PetscInt aoff, i = 0; 352cb1e1211SMatthew G Knepley 353cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(leafSectionAdj, d, &aoff);CHKERRQ(ierr); 354cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 355cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 356cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 357cb1e1211SMatthew G Knepley 358cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 359cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 360cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 361cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 362cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd) { 363cb1e1211SMatthew G Knepley adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd; 364cb1e1211SMatthew G Knepley ++i; 365cb1e1211SMatthew G Knepley } 366cb1e1211SMatthew G Knepley } 367*76185916SToby Isaac for (q = 0; q < anDof; q++) { 368*76185916SToby Isaac adj[aoff+i] = anchorAdj[anOff+q]; 369*76185916SToby Isaac ++i; 370*76185916SToby Isaac } 371cb1e1211SMatthew G Knepley } 372cb1e1211SMatthew G Knepley } 373cb1e1211SMatthew G Knepley /* Debugging */ 374cb1e1211SMatthew G Knepley if (debug) { 375cb1e1211SMatthew G Knepley IS tmp; 376cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Leaf adjacency indices\n");CHKERRQ(ierr); 377cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 378cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 3790a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 380cb1e1211SMatthew G Knepley } 381cb1e1211SMatthew G Knepley /* Gather adjacenct indices to root */ 382cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(rootSectionAdj, &adjSize);CHKERRQ(ierr); 383785e854fSJed Brown ierr = PetscMalloc1(adjSize, &rootAdj);CHKERRQ(ierr); 384cb1e1211SMatthew G Knepley for (r = 0; r < adjSize; ++r) rootAdj[r] = -1; 38547a6104aSMatthew G. Knepley if (doComm) { 386cb1e1211SMatthew G Knepley ierr = PetscSFGatherBegin(sfAdj, MPIU_INT, adj, rootAdj);CHKERRQ(ierr); 387cb1e1211SMatthew G Knepley ierr = PetscSFGatherEnd(sfAdj, MPIU_INT, adj, rootAdj);CHKERRQ(ierr); 388cb1e1211SMatthew G Knepley } 389cb1e1211SMatthew G Knepley ierr = PetscSFDestroy(&sfAdj);CHKERRQ(ierr); 390cb1e1211SMatthew G Knepley ierr = PetscFree(adj);CHKERRQ(ierr); 391cb1e1211SMatthew G Knepley /* Debugging */ 392cb1e1211SMatthew G Knepley if (debug) { 393cb1e1211SMatthew G Knepley IS tmp; 394cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices after gather\n");CHKERRQ(ierr); 395cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 396cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 3970a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 398cb1e1211SMatthew G Knepley } 399cb1e1211SMatthew G Knepley /* Add in local adjacency indices for owned dofs on interface (roots) */ 400cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 401*76185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff; 402cb1e1211SMatthew G Knepley 403cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 404cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 405cb1e1211SMatthew G Knepley if (!dof) continue; 406cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 407cb1e1211SMatthew G Knepley if (adof <= 0) continue; 40870034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 409*76185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 410*76185916SToby Isaac ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr); 411cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 412cb1e1211SMatthew G Knepley PetscInt adof, aoff, i; 413cb1e1211SMatthew G Knepley 414cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr); 415cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr); 416cb1e1211SMatthew G Knepley i = adof-1; 417*76185916SToby Isaac for (q = 0; q < anDof; q++) { 418*76185916SToby Isaac rootAdj[aoff+i] = anchorAdj[anOff+q]; 419*76185916SToby Isaac --i; 420*76185916SToby Isaac } 421cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 422cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 423cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 424cb1e1211SMatthew G Knepley 425cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 426cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 427cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 428cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 429cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd) { 430cb1e1211SMatthew G Knepley rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd; 431cb1e1211SMatthew G Knepley --i; 432cb1e1211SMatthew G Knepley } 433cb1e1211SMatthew G Knepley } 434cb1e1211SMatthew G Knepley } 435cb1e1211SMatthew G Knepley } 436cb1e1211SMatthew G Knepley /* Debugging */ 437cb1e1211SMatthew G Knepley if (debug) { 438cb1e1211SMatthew G Knepley IS tmp; 439cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices\n");CHKERRQ(ierr); 440cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 441cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 4420a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 443cb1e1211SMatthew G Knepley } 444cb1e1211SMatthew G Knepley /* Compress indices */ 445cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr); 446cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 447cb1e1211SMatthew G Knepley PetscInt dof, cdof, off, d; 448cb1e1211SMatthew G Knepley PetscInt adof, aoff; 449cb1e1211SMatthew G Knepley 450cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 451cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 452cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 453cb1e1211SMatthew G Knepley if (!dof) continue; 454cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 455cb1e1211SMatthew G Knepley if (adof <= 0) continue; 456cb1e1211SMatthew G Knepley for (d = off; d < off+dof-cdof; ++d) { 457cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr); 458cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr); 459cb1e1211SMatthew G Knepley ierr = PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]);CHKERRQ(ierr); 460cb1e1211SMatthew G Knepley ierr = PetscSectionSetDof(rootSectionAdj, d, adof);CHKERRQ(ierr); 461cb1e1211SMatthew G Knepley } 462cb1e1211SMatthew G Knepley } 463cb1e1211SMatthew G Knepley /* Debugging */ 464cb1e1211SMatthew G Knepley if (debug) { 465cb1e1211SMatthew G Knepley IS tmp; 466cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n");CHKERRQ(ierr); 467cb1e1211SMatthew G Knepley ierr = PetscSectionView(rootSectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 468cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices after compression\n");CHKERRQ(ierr); 469cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 470cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 4710a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 472cb1e1211SMatthew G Knepley } 473cb1e1211SMatthew G Knepley /* Build adjacency section: Maps global indices to sets of adjacent global indices */ 474cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd);CHKERRQ(ierr); 475cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, §ionAdj);CHKERRQ(ierr); 476cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd);CHKERRQ(ierr); 477cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 478*76185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof; 479cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 480cb1e1211SMatthew G Knepley 481cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 482cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 483cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 484cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 485cb1e1211SMatthew G Knepley for (d = 0; d < dof-cdof; ++d) { 486cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 487cb1e1211SMatthew G Knepley 488cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr); 489cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr); 490cb1e1211SMatthew G Knepley if (ldof > 0) { 491cb1e1211SMatthew G Knepley /* We do not own this point */ 492cb1e1211SMatthew G Knepley } else if (rdof > 0) { 493cb1e1211SMatthew G Knepley ierr = PetscSectionSetDof(sectionAdj, goff+d, rdof);CHKERRQ(ierr); 494cb1e1211SMatthew G Knepley } else { 495cb1e1211SMatthew G Knepley found = PETSC_FALSE; 496cb1e1211SMatthew G Knepley } 497cb1e1211SMatthew G Knepley } 498cb1e1211SMatthew G Knepley if (found) continue; 499cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 500cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 50170034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 502cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 503cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 504cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, noff; 505cb1e1211SMatthew G Knepley 506cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 507cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 508cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 509cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, padj, &noff);CHKERRQ(ierr); 510cb1e1211SMatthew G Knepley for (d = goff; d < goff+dof-cdof; ++d) { 511cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(sectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 512cb1e1211SMatthew G Knepley } 513cb1e1211SMatthew G Knepley } 514*76185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 515*76185916SToby Isaac if (anDof) { 516*76185916SToby Isaac for (d = goff; d < goff+dof-cdof; ++d) { 517*76185916SToby Isaac ierr = PetscSectionAddDof(sectionAdj, d, anDof);CHKERRQ(ierr); 518*76185916SToby Isaac } 519*76185916SToby Isaac } 520cb1e1211SMatthew G Knepley } 521cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(sectionAdj);CHKERRQ(ierr); 522cb1e1211SMatthew G Knepley if (debug) { 523cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency Section for Preallocation:\n");CHKERRQ(ierr); 524cb1e1211SMatthew G Knepley ierr = PetscSectionView(sectionAdj, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 525cb1e1211SMatthew G Knepley } 526cb1e1211SMatthew G Knepley /* Get adjacent indices */ 527cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(sectionAdj, &numCols);CHKERRQ(ierr); 528785e854fSJed Brown ierr = PetscMalloc1(numCols, &cols);CHKERRQ(ierr); 529cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 530*76185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff; 531cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 532cb1e1211SMatthew G Knepley 533cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 534cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 535cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 536cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 537cb1e1211SMatthew G Knepley for (d = 0; d < dof-cdof; ++d) { 538cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 539cb1e1211SMatthew G Knepley 540cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr); 541cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr); 542cb1e1211SMatthew G Knepley if (ldof > 0) { 543cb1e1211SMatthew G Knepley /* We do not own this point */ 544cb1e1211SMatthew G Knepley } else if (rdof > 0) { 545cb1e1211SMatthew G Knepley PetscInt aoff, roff; 546cb1e1211SMatthew G Knepley 547cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, goff+d, &aoff);CHKERRQ(ierr); 548cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, off+d, &roff);CHKERRQ(ierr); 549cb1e1211SMatthew G Knepley ierr = PetscMemcpy(&cols[aoff], &rootAdj[roff], rdof * sizeof(PetscInt));CHKERRQ(ierr); 550cb1e1211SMatthew G Knepley } else { 551cb1e1211SMatthew G Knepley found = PETSC_FALSE; 552cb1e1211SMatthew G Knepley } 553cb1e1211SMatthew G Knepley } 554cb1e1211SMatthew G Knepley if (found) continue; 55570034214SMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &numAdj, &tmpAdj);CHKERRQ(ierr); 556*76185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 557*76185916SToby Isaac ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr); 558cb1e1211SMatthew G Knepley for (d = goff; d < goff+dof-cdof; ++d) { 559cb1e1211SMatthew G Knepley PetscInt adof, aoff, i = 0; 560cb1e1211SMatthew G Knepley 561cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, d, &adof);CHKERRQ(ierr); 562cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, d, &aoff);CHKERRQ(ierr); 563cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 564cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 565cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 566cb1e1211SMatthew G Knepley const PetscInt *ncind; 567cb1e1211SMatthew G Knepley 568cb1e1211SMatthew G Knepley /* Adjacent points may not be in the section chart */ 569cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 570cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 571cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 572cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintIndices(section, padj, &ncind);CHKERRQ(ierr); 573cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 574cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd, ++i) { 575cb1e1211SMatthew G Knepley cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd; 576cb1e1211SMatthew G Knepley } 577cb1e1211SMatthew G Knepley } 578*76185916SToby Isaac for (q = 0; q < anDof; q++, i++) { 579*76185916SToby Isaac cols[aoff+i] = anchorAdj[anOff + q]; 580*76185916SToby Isaac } 581cb1e1211SMatthew 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); 582cb1e1211SMatthew G Knepley } 583cb1e1211SMatthew G Knepley } 584*76185916SToby Isaac ierr = PetscSectionDestroy(&anchorSectionAdj);CHKERRQ(ierr); 585cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(&leafSectionAdj);CHKERRQ(ierr); 586cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(&rootSectionAdj);CHKERRQ(ierr); 587*76185916SToby Isaac ierr = PetscFree(anchorAdj);CHKERRQ(ierr); 588cb1e1211SMatthew G Knepley ierr = PetscFree(rootAdj);CHKERRQ(ierr); 58970034214SMatthew G. Knepley ierr = PetscFree(tmpAdj);CHKERRQ(ierr); 590cb1e1211SMatthew G Knepley /* Debugging */ 591cb1e1211SMatthew G Knepley if (debug) { 592cb1e1211SMatthew G Knepley IS tmp; 593cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Column indices\n");CHKERRQ(ierr); 594cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 595cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 5960a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 597cb1e1211SMatthew G Knepley } 598cb1e1211SMatthew G Knepley /* Create allocation vectors from adjacency graph */ 599cb1e1211SMatthew G Knepley ierr = MatGetLocalSize(A, &locRows, NULL);CHKERRQ(ierr); 600cb1e1211SMatthew G Knepley ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)A), &rLayout);CHKERRQ(ierr); 601cb1e1211SMatthew G Knepley ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 602cb1e1211SMatthew G Knepley ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 603cb1e1211SMatthew G Knepley ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 604cb1e1211SMatthew G Knepley ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 605cb1e1211SMatthew G Knepley ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 606cb1e1211SMatthew G Knepley /* Only loop over blocks of rows */ 607cb1e1211SMatthew 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); 608cb1e1211SMatthew G Knepley for (r = rStart/bs; r < rEnd/bs; ++r) { 609cb1e1211SMatthew G Knepley const PetscInt row = r*bs; 610cb1e1211SMatthew G Knepley PetscInt numCols, cStart, c; 611cb1e1211SMatthew G Knepley 612cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, row, &numCols);CHKERRQ(ierr); 613cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, row, &cStart);CHKERRQ(ierr); 614cb1e1211SMatthew G Knepley for (c = cStart; c < cStart+numCols; ++c) { 615cb1e1211SMatthew G Knepley if ((cols[c] >= rStart*bs) && (cols[c] < rEnd*bs)) { 616cb1e1211SMatthew G Knepley ++dnz[r-rStart]; 617cb1e1211SMatthew G Knepley if (cols[c] >= row) ++dnzu[r-rStart]; 618cb1e1211SMatthew G Knepley } else { 619cb1e1211SMatthew G Knepley ++onz[r-rStart]; 620cb1e1211SMatthew G Knepley if (cols[c] >= row) ++onzu[r-rStart]; 621cb1e1211SMatthew G Knepley } 622cb1e1211SMatthew G Knepley } 623cb1e1211SMatthew G Knepley } 624cb1e1211SMatthew G Knepley if (bs > 1) { 625cb1e1211SMatthew G Knepley for (r = 0; r < locRows/bs; ++r) { 626cb1e1211SMatthew G Knepley dnz[r] /= bs; 627cb1e1211SMatthew G Knepley onz[r] /= bs; 628cb1e1211SMatthew G Knepley dnzu[r] /= bs; 629cb1e1211SMatthew G Knepley onzu[r] /= bs; 630cb1e1211SMatthew G Knepley } 631cb1e1211SMatthew G Knepley } 632cb1e1211SMatthew G Knepley /* Set matrix pattern */ 633cb1e1211SMatthew G Knepley ierr = MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu);CHKERRQ(ierr); 634cb1e1211SMatthew G Knepley ierr = MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 63589545effSMatthew G. Knepley /* Check for symmetric storage */ 63689545effSMatthew G. Knepley ierr = MatGetType(A, &mtype);CHKERRQ(ierr); 63789545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr); 63889545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr); 63989545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr); 64089545effSMatthew G. Knepley if (isSymBlock || isSymSeqBlock || isSymMPIBlock) {ierr = MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);} 641cb1e1211SMatthew G Knepley /* Fill matrix with zeros */ 642cb1e1211SMatthew G Knepley if (fillMatrix) { 643cb1e1211SMatthew G Knepley PetscScalar *values; 644cb1e1211SMatthew G Knepley PetscInt maxRowLen = 0; 645cb1e1211SMatthew G Knepley 646cb1e1211SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 647cb1e1211SMatthew G Knepley PetscInt len; 648cb1e1211SMatthew G Knepley 649cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, r, &len);CHKERRQ(ierr); 650cb1e1211SMatthew G Knepley maxRowLen = PetscMax(maxRowLen, len); 651cb1e1211SMatthew G Knepley } 6521795a4d1SJed Brown ierr = PetscCalloc1(maxRowLen, &values);CHKERRQ(ierr); 653cb1e1211SMatthew G Knepley for (r = rStart; r < rEnd; ++r) { 654cb1e1211SMatthew G Knepley PetscInt numCols, cStart; 655cb1e1211SMatthew G Knepley 656cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr); 657cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr); 658cb1e1211SMatthew G Knepley ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr); 659cb1e1211SMatthew G Knepley } 660cb1e1211SMatthew G Knepley ierr = PetscFree(values);CHKERRQ(ierr); 661cb1e1211SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 662cb1e1211SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 663cb1e1211SMatthew G Knepley } 6648b0b4c70SToby Isaac /* restore original useConstraints */ 6658b0b4c70SToby Isaac ierr = DMPlexSetAdjacencyUseConstraints(dm,useConstraints);CHKERRQ(ierr); 666cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(§ionAdj);CHKERRQ(ierr); 667cb1e1211SMatthew G Knepley ierr = PetscFree(cols);CHKERRQ(ierr); 668a72f3261SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr); 669cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 670cb1e1211SMatthew G Knepley } 671cb1e1211SMatthew G Knepley 672cb1e1211SMatthew G Knepley #if 0 673cb1e1211SMatthew G Knepley #undef __FUNCT__ 674cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexPreallocateOperator_2" 675cb1e1211SMatthew 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) 676cb1e1211SMatthew G Knepley { 677cb1e1211SMatthew G Knepley PetscInt *tmpClosure,*tmpAdj,*visits; 678cb1e1211SMatthew G Knepley PetscInt c,cStart,cEnd,pStart,pEnd; 679cb1e1211SMatthew G Knepley PetscErrorCode ierr; 680cb1e1211SMatthew G Knepley 681cb1e1211SMatthew G Knepley PetscFunctionBegin; 682c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 683cb1e1211SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 684cb1e1211SMatthew G Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 685cb1e1211SMatthew G Knepley 686e7b80042SMatthew G. Knepley maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1)); 687cb1e1211SMatthew G Knepley 688cb1e1211SMatthew G Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 689cb1e1211SMatthew G Knepley npoints = pEnd - pStart; 690cb1e1211SMatthew G Knepley 691dcca6d9dSJed Brown ierr = PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits);CHKERRQ(ierr); 692cb1e1211SMatthew G Knepley ierr = PetscMemzero(lvisits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr); 693cb1e1211SMatthew G Knepley ierr = PetscMemzero(visits,(pEnd-pStart)*sizeof(PetscInt));CHKERRQ(ierr); 694cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 695cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 696cb1e1211SMatthew G Knepley PetscInt *support = tmpClosure; 697cb1e1211SMatthew G Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support);CHKERRQ(ierr); 698cb1e1211SMatthew G Knepley for (p=0; p<supportSize; p++) lvisits[support[p]]++; 699cb1e1211SMatthew G Knepley } 700cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr); 701cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd (sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr); 702cb1e1211SMatthew G Knepley ierr = PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr); 703cb1e1211SMatthew G Knepley ierr = PetscSFBcastEnd (sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr); 704cb1e1211SMatthew G Knepley 705cb1e1211SMatthew G Knepley ierr = PetscSFGetRanks();CHKERRQ(ierr); 706cb1e1211SMatthew G Knepley 707cb1e1211SMatthew G Knepley 708dcca6d9dSJed Brown ierr = PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner);CHKERRQ(ierr); 709cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 710cb1e1211SMatthew G Knepley ierr = PetscMemzero(cellmat,maxClosureSize*maxClosureSize*sizeof(PetscInt));CHKERRQ(ierr); 711cb1e1211SMatthew G Knepley /* 712cb1e1211SMatthew G Knepley Depth-first walk of transitive closure. 713cb1e1211SMatthew 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. 714cb1e1211SMatthew G Knepley This contribution is added to dnz if owning ranks of p and q match, to onz otherwise. 715cb1e1211SMatthew G Knepley */ 716cb1e1211SMatthew G Knepley } 717cb1e1211SMatthew G Knepley 718cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM);CHKERRQ(ierr); 719cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd (sf,MPIU_INT,lonz,onz,MPI_SUM);CHKERRQ(ierr); 720cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 721cb1e1211SMatthew G Knepley } 722cb1e1211SMatthew G Knepley #endif 723