1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/isimpl.h> 30c312b8eSJed Brown #include <petscsf.h> 4a9fb6443SMatthew G. Knepley #include <petscds.h> 5cb1e1211SMatthew G Knepley 676185916SToby Isaac /* get adjacencies due to point-to-point constraints that can't be found with DMPlexGetAdjacency() */ 7d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[]) 8d71ae5a4SJacob Faibussowitsch { 976185916SToby Isaac PetscInt pStart, pEnd; 10e101f074SMatthew G. Knepley PetscSection section, sectionGlobal, adjSec, aSec; 1176185916SToby Isaac IS aIS; 1276185916SToby Isaac 1376185916SToby Isaac PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 159566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, §ionGlobal)); 169566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)section), &adjSec)); 179566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 189566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(adjSec, pStart, pEnd)); 1976185916SToby Isaac 209566063dSJacob Faibussowitsch PetscCall(DMPlexGetAnchors(dm, &aSec, &aIS)); 2176185916SToby Isaac if (aSec) { 2276185916SToby Isaac const PetscInt *anchors; 2376185916SToby Isaac PetscInt p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize; 2476185916SToby Isaac PetscInt *tmpAdjP = NULL, *tmpAdjQ = NULL; 2576185916SToby Isaac PetscSection inverseSec; 2676185916SToby Isaac 2776185916SToby Isaac /* invert the constraint-to-anchor map */ 289566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)aSec), &inverseSec)); 299566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(inverseSec, pStart, pEnd)); 309566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(aIS, &aSize)); 319566063dSJacob Faibussowitsch PetscCall(ISGetIndices(aIS, &anchors)); 3276185916SToby Isaac 3376185916SToby Isaac for (p = 0; p < aSize; p++) { 3476185916SToby Isaac PetscInt a = anchors[p]; 3576185916SToby Isaac 369566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(inverseSec, a, 1)); 3776185916SToby Isaac } 389566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(inverseSec)); 399566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(inverseSec, &iSize)); 409566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(iSize, &inverse)); 419566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(pEnd - pStart, &offsets)); 429566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(aSec, &aStart, &aEnd)); 4376185916SToby Isaac for (p = aStart; p < aEnd; p++) { 4476185916SToby Isaac PetscInt dof, off; 4576185916SToby Isaac 469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(aSec, p, &dof)); 479566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(aSec, p, &off)); 4876185916SToby Isaac 4976185916SToby Isaac for (q = 0; q < dof; q++) { 5076185916SToby Isaac PetscInt iOff; 5176185916SToby Isaac 5276185916SToby Isaac a = anchors[off + q]; 539566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(inverseSec, a, &iOff)); 5476185916SToby Isaac inverse[iOff + offsets[a - pStart]++] = p; 5576185916SToby Isaac } 5676185916SToby Isaac } 579566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(aIS, &anchors)); 589566063dSJacob Faibussowitsch PetscCall(PetscFree(offsets)); 5976185916SToby Isaac 6076185916SToby Isaac /* construct anchorAdj and adjSec 6176185916SToby Isaac * 6276185916SToby Isaac * loop over anchors: 6376185916SToby Isaac * construct anchor adjacency 6476185916SToby Isaac * loop over constrained: 6576185916SToby Isaac * construct constrained adjacency 6676185916SToby Isaac * if not in anchor adjacency, add to dofs 6776185916SToby Isaac * setup adjSec, allocate anchorAdj 6876185916SToby Isaac * loop over anchors: 6976185916SToby Isaac * construct anchor adjacency 7076185916SToby Isaac * loop over constrained: 7176185916SToby Isaac * construct constrained adjacency 7276185916SToby Isaac * if not in anchor adjacency 7376185916SToby Isaac * if not already in list, put in list 7476185916SToby Isaac * sort, unique, reduce dof count 7576185916SToby Isaac * optional: compactify 7676185916SToby Isaac */ 7776185916SToby Isaac for (p = pStart; p < pEnd; p++) { 7876185916SToby Isaac PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE; 7976185916SToby Isaac 809566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(inverseSec, p, &iDof)); 8176185916SToby Isaac if (!iDof) continue; 829566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff)); 839566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP)); 8476185916SToby Isaac for (i = 0; i < iDof; i++) { 8576185916SToby Isaac PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE; 8676185916SToby Isaac 8776185916SToby Isaac q = inverse[iOff + i]; 889566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ)); 8976185916SToby Isaac for (r = 0; r < numAdjQ; r++) { 9076185916SToby Isaac qAdj = tmpAdjQ[r]; 9176185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue; 9276185916SToby Isaac for (s = 0; s < numAdjP; s++) { 9376185916SToby Isaac if (qAdj == tmpAdjP[s]) break; 9476185916SToby Isaac } 9576185916SToby Isaac if (s < numAdjP) continue; 969566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof)); 979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof)); 9876185916SToby Isaac iNew += qAdjDof - qAdjCDof; 9976185916SToby Isaac } 1009566063dSJacob Faibussowitsch PetscCall(PetscSectionAddDof(adjSec, p, iNew)); 10176185916SToby Isaac } 10276185916SToby Isaac } 10376185916SToby Isaac 1049566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(adjSec)); 1059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(adjSec, &adjSize)); 1069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(adjSize, &adj)); 10776185916SToby Isaac 10876185916SToby Isaac for (p = pStart; p < pEnd; p++) { 10976185916SToby Isaac PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE; 11076185916SToby Isaac 1119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(inverseSec, p, &iDof)); 11276185916SToby Isaac if (!iDof) continue; 1139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(inverseSec, p, &iOff)); 1149566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, PETSC_TRUE, &numAdjP, &tmpAdjP)); 1159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(adjSec, p, &aDof)); 1169566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(adjSec, p, &aOff)); 11776185916SToby Isaac aOffOrig = aOff; 11876185916SToby Isaac for (i = 0; i < iDof; i++) { 11976185916SToby Isaac PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE; 12076185916SToby Isaac 12176185916SToby Isaac q = inverse[iOff + i]; 1229566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, q, useCone, useClosure, PETSC_TRUE, &numAdjQ, &tmpAdjQ)); 12376185916SToby Isaac for (r = 0; r < numAdjQ; r++) { 12476185916SToby Isaac qAdj = tmpAdjQ[r]; 12576185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue; 12676185916SToby Isaac for (s = 0; s < numAdjP; s++) { 12776185916SToby Isaac if (qAdj == tmpAdjP[s]) break; 12876185916SToby Isaac } 12976185916SToby Isaac if (s < numAdjP) continue; 1309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, qAdj, &qAdjDof)); 1319566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, qAdj, &qAdjCDof)); 1329566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, qAdj, &qAdjOff)); 133ad540459SPierre Jolivet for (nd = 0; nd < qAdjDof - qAdjCDof; ++nd) adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff + 1) : qAdjOff) + nd; 13476185916SToby Isaac } 13576185916SToby Isaac } 1368e3a54c0SPierre Jolivet PetscCall(PetscSortRemoveDupsInt(&aDof, PetscSafePointerPlusOffset(adj, aOffOrig))); 1379566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(adjSec, p, aDof)); 13876185916SToby Isaac } 13976185916SToby Isaac *anchorAdj = adj; 14076185916SToby Isaac 14176185916SToby Isaac /* clean up */ 1429566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&inverseSec)); 1439566063dSJacob Faibussowitsch PetscCall(PetscFree(inverse)); 1449566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpAdjP)); 1459566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpAdjQ)); 1469371c9d4SSatish Balay } else { 14776185916SToby Isaac *anchorAdj = NULL; 1489566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(adjSec)); 14976185916SToby Isaac } 15076185916SToby Isaac *anchorSectionAdj = adjSec; 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15276185916SToby Isaac } 15376185916SToby Isaac 154*f0b52e9eSJames Wright // Based off of `PetscIntViewNumColumns()` 155*f0b52e9eSJames Wright static PetscErrorCode PetscIntViewPairs(PetscInt N, PetscInt Ncol, const PetscInt idx1[], const PetscInt idx2[], PetscViewer viewer) 156*f0b52e9eSJames Wright { 157*f0b52e9eSJames Wright PetscMPIInt rank, size; 158*f0b52e9eSJames Wright PetscInt j, i, n = N / Ncol, p = N % Ncol; 159*f0b52e9eSJames Wright PetscBool isascii; 160*f0b52e9eSJames Wright MPI_Comm comm; 161*f0b52e9eSJames Wright 162*f0b52e9eSJames Wright PetscFunctionBegin; 163*f0b52e9eSJames Wright if (!viewer) viewer = PETSC_VIEWER_STDOUT_SELF; 164*f0b52e9eSJames Wright if (N) PetscAssertPointer(idx1, 3); 165*f0b52e9eSJames Wright if (N) PetscAssertPointer(idx2, 4); 166*f0b52e9eSJames Wright PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 5); 167*f0b52e9eSJames Wright PetscCall(PetscObjectGetComm((PetscObject)viewer, &comm)); 168*f0b52e9eSJames Wright PetscCallMPI(MPI_Comm_size(comm, &size)); 169*f0b52e9eSJames Wright PetscCallMPI(MPI_Comm_rank(comm, &rank)); 170*f0b52e9eSJames Wright 171*f0b52e9eSJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 172*f0b52e9eSJames Wright if (isascii) { 173*f0b52e9eSJames Wright PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 174*f0b52e9eSJames Wright for (i = 0; i < n; i++) { 175*f0b52e9eSJames Wright if (size > 1) { 176*f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT ":", rank, Ncol * i)); 177*f0b52e9eSJames Wright } else { 178*f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT ":", Ncol * i)); 179*f0b52e9eSJames Wright } 180*f0b52e9eSJames Wright for (j = 0; j < Ncol; j++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%" PetscInt_FMT ", %" PetscInt_FMT ")", idx1[i * Ncol + j], idx2[i * Ncol + j])); 181*f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 182*f0b52e9eSJames Wright } 183*f0b52e9eSJames Wright if (p) { 184*f0b52e9eSJames Wright if (size > 1) { 185*f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "[%d] %" PetscInt_FMT ":", rank, Ncol * n)); 186*f0b52e9eSJames Wright } else { 187*f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%" PetscInt_FMT ":", Ncol * n)); 188*f0b52e9eSJames Wright } 189*f0b52e9eSJames Wright for (i = 0; i < p; i++) PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%" PetscInt_FMT ", %" PetscInt_FMT ")", idx1[Ncol * n + i], idx2[Ncol * n + i])); 190*f0b52e9eSJames Wright PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 191*f0b52e9eSJames Wright } 192*f0b52e9eSJames Wright PetscCall(PetscViewerFlush(viewer)); 193*f0b52e9eSJames Wright PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 194*f0b52e9eSJames Wright } else { 195*f0b52e9eSJames Wright const char *tname; 196*f0b52e9eSJames Wright PetscCall(PetscObjectGetName((PetscObject)viewer, &tname)); 197*f0b52e9eSJames Wright SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle that PetscViewer of type %s", tname); 198*f0b52e9eSJames Wright } 199*f0b52e9eSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 200*f0b52e9eSJames Wright } 201*f0b52e9eSJames Wright 2028c5add6aSPierre Jolivet // Determine if any of the local adjacencies match a leaf and root of the pointSF. 203a38eeca9SJames Wright // When using isoperiodic boundary conditions, it is possible for periodic (leaf) and donor (root) pairs to be on the same rank. 2048c5add6aSPierre Jolivet // This check is done to ensure the adjacency in these cases is only counted for one of the mesh points rather than both. 205a38eeca9SJames Wright static inline PetscErrorCode AdjancencyContainsLeafRootPair(PetscInt num_pairs, const PetscInt leaves[], const PetscInt roots[], PetscInt numAdj, const PetscInt tmpAdj[], PetscInt *num_leaves_found, PetscInt leaves_found[]) 206a38eeca9SJames Wright { 207a38eeca9SJames Wright PetscInt root_index = -1, leaf_, num_roots_found = 0; 208a38eeca9SJames Wright 209a38eeca9SJames Wright PetscFunctionBeginUser; 210a38eeca9SJames Wright *num_leaves_found = 0; 211a38eeca9SJames Wright for (PetscInt q = 0; q < numAdj; q++) { 212a38eeca9SJames Wright const PetscInt padj = tmpAdj[q]; 213a38eeca9SJames Wright PetscCall(PetscFindInt(padj, num_pairs, roots, &root_index)); 214a38eeca9SJames Wright if (root_index >= 0) { 215a38eeca9SJames Wright leaves_found[num_roots_found] = root_index; // Initially use leaves_found to store pair indices 216a38eeca9SJames Wright num_roots_found++; 217a38eeca9SJames Wright break; 218a38eeca9SJames Wright } 219a38eeca9SJames Wright } 220a38eeca9SJames Wright if (num_roots_found == 0) PetscFunctionReturn(PETSC_SUCCESS); 221a38eeca9SJames Wright 222a38eeca9SJames Wright for (PetscInt i = 0; i < num_roots_found; i++) { 223a38eeca9SJames Wright leaf_ = leaves[root_index]; 224a38eeca9SJames Wright for (PetscInt q = 0; q < numAdj; q++) { 225a38eeca9SJames Wright if (tmpAdj[q] == leaf_) { 226a38eeca9SJames Wright leaves_found[*num_leaves_found] = leaf_; 227a38eeca9SJames Wright (*num_leaves_found)++; 228a38eeca9SJames Wright continue; 229a38eeca9SJames Wright } 230a38eeca9SJames Wright } 231a38eeca9SJames Wright } 232a38eeca9SJames Wright 233a38eeca9SJames Wright PetscCall(PetscIntSortSemiOrdered(*num_leaves_found, leaves_found)); 234a38eeca9SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 235a38eeca9SJames Wright } 236a38eeca9SJames Wright 237d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx) 238d71ae5a4SJacob Faibussowitsch { 239cb1e1211SMatthew G Knepley MPI_Comm comm; 240a38eeca9SJames Wright PetscMPIInt myrank; 241a9fb6443SMatthew G. Knepley PetscBool doCommLocal, doComm, debug = PETSC_FALSE; 242a9fb6443SMatthew G. Knepley PetscSF sf, sfAdj; 243e101f074SMatthew G. Knepley PetscSection section, sectionGlobal, leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj; 244a9fb6443SMatthew G. Knepley PetscInt nroots, nleaves, l, p, r; 245cb1e1211SMatthew G Knepley const PetscInt *leaves; 246cb1e1211SMatthew G Knepley const PetscSFNode *remotes; 247cb1e1211SMatthew G Knepley PetscInt dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols; 248a38eeca9SJames Wright PetscInt *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets, *rootsMyRankPair = NULL, *leavesMyRankPair = NULL; 249a38eeca9SJames Wright PetscInt adjSize, numMyRankPair = 0; 250cb1e1211SMatthew G Knepley 251cb1e1211SMatthew G Knepley PetscFunctionBegin; 2529566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 2539566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL)); 254a38eeca9SJames Wright PetscCallMPI(MPI_Comm_rank(comm, &myrank)); 2559566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 256a38eeca9SJames Wright PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf)); 2579566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 2589566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, §ionGlobal)); 2599566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 260a38eeca9SJames Wright doCommLocal = nroots >= 0 ? PETSC_TRUE : PETSC_FALSE; 2615440e5dcSBarry Smith PetscCallMPI(MPIU_Allreduce(&doCommLocal, &doComm, 1, MPI_C_BOOL, MPI_LAND, comm)); 262cb1e1211SMatthew G Knepley /* Create section for dof adjacency (dof ==> # adj dof) */ 2639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 2649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(section, &numDof)); 2659566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &leafSectionAdj)); 2669566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(leafSectionAdj, 0, numDof)); 2679566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, &rootSectionAdj)); 2689566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(rootSectionAdj, 0, numDof)); 2699566063dSJacob Faibussowitsch PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes)); 270a38eeca9SJames Wright 271a38eeca9SJames Wright // Store leaf-root pairs if remote.rank is the current rank 272a38eeca9SJames Wright if (nleaves >= 0) PetscCall(PetscMalloc2(nleaves, &rootsMyRankPair, nleaves, &leavesMyRankPair)); 273a38eeca9SJames Wright for (PetscInt l = 0; l < nleaves; l++) { 274a38eeca9SJames Wright if (remotes[l].rank == myrank) { 275a38eeca9SJames Wright rootsMyRankPair[numMyRankPair] = remotes[l].index; 276a38eeca9SJames Wright leavesMyRankPair[numMyRankPair] = leaves[l]; 277a38eeca9SJames Wright numMyRankPair++; 278a38eeca9SJames Wright } 279a38eeca9SJames Wright } 280a38eeca9SJames Wright PetscCall(PetscIntSortSemiOrdered(numMyRankPair, rootsMyRankPair)); 281a38eeca9SJames Wright PetscCall(PetscIntSortSemiOrdered(numMyRankPair, leavesMyRankPair)); 282a38eeca9SJames Wright if (debug) { 283*f0b52e9eSJames Wright PetscCall(PetscPrintf(comm, "Root/leaf pairs on the same rank:\n")); 284*f0b52e9eSJames Wright PetscCall(PetscIntViewPairs(numMyRankPair, 5, rootsMyRankPair, leavesMyRankPair, NULL)); 285a38eeca9SJames Wright } 286cb1e1211SMatthew G Knepley /* 287964bf7afSMatthew G. Knepley section - maps points to (# dofs, local dofs) 288964bf7afSMatthew G. Knepley sectionGlobal - maps points to (# dofs, global dofs) 289964bf7afSMatthew G. Knepley leafSectionAdj - maps unowned local dofs to # adj dofs 290964bf7afSMatthew G. Knepley rootSectionAdj - maps owned local dofs to # adj dofs 291964bf7afSMatthew G. Knepley adj - adj global dofs indexed by leafSectionAdj 292964bf7afSMatthew G. Knepley rootAdj - adj global dofs indexed by rootSectionAdj 293964bf7afSMatthew G. Knepley sf - describes shared points across procs 294964bf7afSMatthew G. Knepley sfDof - describes shared dofs across procs 295964bf7afSMatthew G. Knepley sfAdj - describes shared adjacent dofs across procs 296cb1e1211SMatthew G Knepley ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point. 29776185916SToby Isaac (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj 29876185916SToby Isaac (This is done in DMPlexComputeAnchorAdjacencies()) 299cb1e1211SMatthew G Knepley 1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj 300cb1e1211SMatthew G Knepley Reduce those counts to rootSectionAdj (now redundantly counting some interface points) 301cb1e1211SMatthew G Knepley 2. Visit owned points on interface, count adjacencies placing in rootSectionAdj 302cb1e1211SMatthew G Knepley Create sfAdj connecting rootSectionAdj and leafSectionAdj 303cb1e1211SMatthew G Knepley 3. Visit unowned points on interface, write adjacencies to adj 304cb1e1211SMatthew G Knepley Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies) 305cb1e1211SMatthew G Knepley 4. Visit owned points on interface, write adjacencies to rootAdj 306cb1e1211SMatthew G Knepley Remove redundancy in rootAdj 307cb1e1211SMatthew G Knepley ** The last two traversals use transitive closure 308cb1e1211SMatthew G Knepley 5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj) 309cb1e1211SMatthew G Knepley Allocate memory addressed by sectionAdj (cols) 310cb1e1211SMatthew G Knepley 6. Visit all owned points in the subdomain, insert dof adjacencies into cols 311cb1e1211SMatthew G Knepley ** Knowing all the column adjacencies, check ownership and sum into dnz and onz 312cb1e1211SMatthew G Knepley */ 3139566063dSJacob Faibussowitsch PetscCall(DMPlexComputeAnchorAdjacencies(dm, useCone, useClosure, &anchorSectionAdj, &anchorAdj)); 314cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 31576185916SToby Isaac PetscInt dof, off, d, q, anDof; 31670034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 317cb1e1211SMatthew G Knepley 318fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 3199566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 3209566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 3219566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj)); 322cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 323cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 324cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 325cb1e1211SMatthew G Knepley 326cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 3279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof)); 3289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof)); 32948a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, ndof - ncdof)); 330cb1e1211SMatthew G Knepley } 3319566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof)); 33276185916SToby Isaac if (anDof) { 33348a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(leafSectionAdj, d, anDof)); 33476185916SToby Isaac } 335cb1e1211SMatthew G Knepley } 3369566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSectionAdj)); 337cb1e1211SMatthew G Knepley if (debug) { 3389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n")); 3399566063dSJacob Faibussowitsch PetscCall(PetscSectionView(leafSectionAdj, NULL)); 340cb1e1211SMatthew G Knepley } 341cb1e1211SMatthew G Knepley /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */ 34247a6104aSMatthew G. Knepley if (doComm) { 3439566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM)); 3449566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM)); 34569c11d05SVaclav Hapla PetscCall(PetscSectionInvalidateMaxDof_Internal(rootSectionAdj)); 346cb1e1211SMatthew G Knepley } 347cb1e1211SMatthew G Knepley if (debug) { 348145b44c9SPierre Jolivet PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots:\n")); 3499566063dSJacob Faibussowitsch PetscCall(PetscSectionView(rootSectionAdj, NULL)); 350cb1e1211SMatthew G Knepley } 351cb1e1211SMatthew G Knepley /* Add in local adjacency sizes for owned dofs on interface (roots) */ 352cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 35376185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof; 354cb1e1211SMatthew G Knepley 3559566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 3569566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 357cb1e1211SMatthew G Knepley if (!dof) continue; 3589566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof)); 359cb1e1211SMatthew G Knepley if (adof <= 0) continue; 3609566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj)); 361cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 362cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 363cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 364cb1e1211SMatthew G Knepley 365cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 3669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof)); 3679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof)); 36848a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, ndof - ncdof)); 369cb1e1211SMatthew G Knepley } 3709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof)); 37176185916SToby Isaac if (anDof) { 37248a46eb9SPierre Jolivet for (d = off; d < off + dof; ++d) PetscCall(PetscSectionAddDof(rootSectionAdj, d, anDof)); 37376185916SToby Isaac } 374cb1e1211SMatthew G Knepley } 3759566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(rootSectionAdj)); 376cb1e1211SMatthew G Knepley if (debug) { 377145b44c9SPierre Jolivet PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots after local additions:\n")); 3789566063dSJacob Faibussowitsch PetscCall(PetscSectionView(rootSectionAdj, NULL)); 379cb1e1211SMatthew G Knepley } 380cb1e1211SMatthew G Knepley /* Create adj SF based on dof SF */ 3819566063dSJacob Faibussowitsch PetscCall(PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets)); 3829566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj)); 3839566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets)); 3846a1a2f7bSJames Wright if (debug) { 3859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Adjacency SF for Preallocation:\n")); 3869566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfAdj, NULL)); 387cb1e1211SMatthew G Knepley } 388cb1e1211SMatthew G Knepley /* Create leaf adjacency */ 3899566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(leafSectionAdj)); 3909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(leafSectionAdj, &adjSize)); 3919566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(adjSize, &adj)); 392cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 39376185916SToby Isaac PetscInt dof, off, d, q, anDof, anOff; 39470034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 395cb1e1211SMatthew G Knepley 396fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 3979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 3989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 3999566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj)); 4009566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof)); 4019566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff)); 402cb1e1211SMatthew G Knepley for (d = off; d < off + dof; ++d) { 403cb1e1211SMatthew G Knepley PetscInt aoff, i = 0; 404cb1e1211SMatthew G Knepley 4059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(leafSectionAdj, d, &aoff)); 406cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 407cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 408cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 409cb1e1211SMatthew G Knepley 410cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 4119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof)); 4129566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof)); 4139566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff)); 414cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof - ncdof; ++nd) { 415cb1e1211SMatthew G Knepley adj[aoff + i] = (ngoff < 0 ? -(ngoff + 1) : ngoff) + nd; 416cb1e1211SMatthew G Knepley ++i; 417cb1e1211SMatthew G Knepley } 418cb1e1211SMatthew G Knepley } 41976185916SToby Isaac for (q = 0; q < anDof; q++) { 42076185916SToby Isaac adj[aoff + i] = anchorAdj[anOff + q]; 42176185916SToby Isaac ++i; 42276185916SToby Isaac } 423cb1e1211SMatthew G Knepley } 424cb1e1211SMatthew G Knepley } 425cb1e1211SMatthew G Knepley /* Debugging */ 426cb1e1211SMatthew G Knepley if (debug) { 4279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Leaf adjacency indices\n")); 4286a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(leafSectionAdj, adj, PETSC_INT, NULL)); 429cb1e1211SMatthew G Knepley } 430543482b8SMatthew G. Knepley /* Gather adjacent indices to root */ 4319566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(rootSectionAdj, &adjSize)); 4329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(adjSize, &rootAdj)); 433cb1e1211SMatthew G Knepley for (r = 0; r < adjSize; ++r) rootAdj[r] = -1; 43447a6104aSMatthew G. Knepley if (doComm) { 435543482b8SMatthew G. Knepley const PetscInt *indegree; 436543482b8SMatthew G. Knepley PetscInt *remoteadj, radjsize = 0; 437543482b8SMatthew G. Knepley 4389566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeBegin(sfAdj, &indegree)); 4399566063dSJacob Faibussowitsch PetscCall(PetscSFComputeDegreeEnd(sfAdj, &indegree)); 440543482b8SMatthew G. Knepley for (p = 0; p < adjSize; ++p) radjsize += indegree[p]; 4419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(radjsize, &remoteadj)); 4429566063dSJacob Faibussowitsch PetscCall(PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj)); 4439566063dSJacob Faibussowitsch PetscCall(PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj)); 4446dba2905SMatthew G. Knepley for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p - 1])) { 445543482b8SMatthew G. Knepley PetscInt s; 4466dba2905SMatthew G. Knepley for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l + s] = remoteadj[r]; 447543482b8SMatthew G. Knepley } 44863a3b9bcSJacob Faibussowitsch PetscCheck(r == radjsize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, r, radjsize); 44963a3b9bcSJacob Faibussowitsch PetscCheck(l == adjSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %" PetscInt_FMT " != %" PetscInt_FMT, l, adjSize); 4509566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteadj)); 451cb1e1211SMatthew G Knepley } 4529566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfAdj)); 4539566063dSJacob Faibussowitsch PetscCall(PetscFree(adj)); 454cb1e1211SMatthew G Knepley /* Debugging */ 455cb1e1211SMatthew G Knepley if (debug) { 4569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Root adjacency indices after gather\n")); 4576a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(rootSectionAdj, rootAdj, PETSC_INT, NULL)); 458cb1e1211SMatthew G Knepley } 459cb1e1211SMatthew G Knepley /* Add in local adjacency indices for owned dofs on interface (roots) */ 460cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 46176185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff; 462cb1e1211SMatthew G Knepley 4639566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 4649566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 465cb1e1211SMatthew G Knepley if (!dof) continue; 4669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof)); 467cb1e1211SMatthew G Knepley if (adof <= 0) continue; 4689566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj)); 4699566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof)); 4709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff)); 471cb1e1211SMatthew G Knepley for (d = off; d < off + dof; ++d) { 472cb1e1211SMatthew G Knepley PetscInt adof, aoff, i; 473cb1e1211SMatthew G Knepley 4749566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof)); 4759566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff)); 476cb1e1211SMatthew G Knepley i = adof - 1; 47776185916SToby Isaac for (q = 0; q < anDof; q++) { 47876185916SToby Isaac rootAdj[aoff + i] = anchorAdj[anOff + q]; 47976185916SToby Isaac --i; 48076185916SToby Isaac } 481cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 482cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 483cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 484cb1e1211SMatthew G Knepley 485cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 4869566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof)); 4879566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof)); 4889566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff)); 489cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof - ncdof; ++nd) { 490cb1e1211SMatthew G Knepley rootAdj[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd; 491cb1e1211SMatthew G Knepley --i; 492cb1e1211SMatthew G Knepley } 493cb1e1211SMatthew G Knepley } 494cb1e1211SMatthew G Knepley } 495cb1e1211SMatthew G Knepley } 496cb1e1211SMatthew G Knepley /* Debugging */ 497cb1e1211SMatthew G Knepley if (debug) { 4986a1a2f7bSJames Wright PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots before compression:\n")); 4996a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(rootSectionAdj, rootAdj, PETSC_INT, NULL)); 500cb1e1211SMatthew G Knepley } 501cb1e1211SMatthew G Knepley /* Compress indices */ 5029566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(rootSectionAdj)); 503cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 504cb1e1211SMatthew G Knepley PetscInt dof, cdof, off, d; 505cb1e1211SMatthew G Knepley PetscInt adof, aoff; 506cb1e1211SMatthew G Knepley 5079566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 5089566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 5099566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 510cb1e1211SMatthew G Knepley if (!dof) continue; 5119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off, &adof)); 512cb1e1211SMatthew G Knepley if (adof <= 0) continue; 513cb1e1211SMatthew G Knepley for (d = off; d < off + dof - cdof; ++d) { 5149566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, d, &adof)); 5159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(rootSectionAdj, d, &aoff)); 5169566063dSJacob Faibussowitsch PetscCall(PetscSortRemoveDupsInt(&adof, &rootAdj[aoff])); 5179566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(rootSectionAdj, d, adof)); 518cb1e1211SMatthew G Knepley } 519cb1e1211SMatthew G Knepley } 520cb1e1211SMatthew G Knepley /* Debugging */ 521cb1e1211SMatthew G Knepley if (debug) { 522145b44c9SPierre Jolivet PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation on Roots after compression:\n")); 5236a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(rootSectionAdj, rootAdj, PETSC_INT, NULL)); 524cb1e1211SMatthew G Knepley } 525cb1e1211SMatthew G Knepley /* Build adjacency section: Maps global indices to sets of adjacent global indices */ 5269566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd)); 5279566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(comm, §ionAdj)); 5289566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd)); 529a38eeca9SJames Wright 530a38eeca9SJames Wright PetscInt *exclude_leaves, num_exclude_leaves, max_adjacency_size = 0; 531a38eeca9SJames Wright PetscCall(DMPlexGetMaxAdjacencySize_Internal(dm, useAnchors, &max_adjacency_size)); 532a38eeca9SJames Wright PetscCall(PetscMalloc1(max_adjacency_size, &exclude_leaves)); 533a38eeca9SJames Wright 534cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 53576185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof; 536cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 537cb1e1211SMatthew G Knepley 5389566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 5399566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 5409566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 5419566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff)); 542cb1e1211SMatthew G Knepley for (d = 0; d < dof - cdof; ++d) { 543cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 544cb1e1211SMatthew G Knepley 5459566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof)); 5469566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof)); 547cb1e1211SMatthew G Knepley if (ldof > 0) { 548cb1e1211SMatthew G Knepley /* We do not own this point */ 549cb1e1211SMatthew G Knepley } else if (rdof > 0) { 5509566063dSJacob Faibussowitsch PetscCall(PetscSectionSetDof(sectionAdj, goff + d, rdof)); 551cb1e1211SMatthew G Knepley } else { 552cb1e1211SMatthew G Knepley found = PETSC_FALSE; 553cb1e1211SMatthew G Knepley } 554cb1e1211SMatthew G Knepley } 555cb1e1211SMatthew G Knepley if (found) continue; 5569566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 5579566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff)); 5589566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj)); 559a38eeca9SJames Wright PetscCall(AdjancencyContainsLeafRootPair(numMyRankPair, leavesMyRankPair, rootsMyRankPair, numAdj, tmpAdj, &num_exclude_leaves, exclude_leaves)); 560cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 561cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 562a38eeca9SJames Wright PetscInt ndof, ncdof, noff, count; 563cb1e1211SMatthew G Knepley 564cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 5659566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof)); 5669566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof)); 5679566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, padj, &noff)); 568a38eeca9SJames Wright // If leaf-root pair are both on this rank, only count root 569a38eeca9SJames Wright PetscCall(PetscFindInt(padj, num_exclude_leaves, exclude_leaves, &count)); 570a38eeca9SJames Wright if (count >= 0) continue; 57148a46eb9SPierre Jolivet for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, ndof - ncdof)); 572cb1e1211SMatthew G Knepley } 5739566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof)); 57476185916SToby Isaac if (anDof) { 57548a46eb9SPierre Jolivet for (d = goff; d < goff + dof - cdof; ++d) PetscCall(PetscSectionAddDof(sectionAdj, d, anDof)); 57676185916SToby Isaac } 577cb1e1211SMatthew G Knepley } 5789566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(sectionAdj)); 579cb1e1211SMatthew G Knepley if (debug) { 5809566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Adjacency Section for Preallocation:\n")); 5819566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sectionAdj, NULL)); 582cb1e1211SMatthew G Knepley } 583cb1e1211SMatthew G Knepley /* Get adjacent indices */ 5849566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(sectionAdj, &numCols)); 5859566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(numCols, &cols)); 586cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 58776185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff; 588cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 589cb1e1211SMatthew G Knepley 5909566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, p, &dof)); 5919566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 5929566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(section, p, &off)); 5939566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, p, &goff)); 594cb1e1211SMatthew G Knepley for (d = 0; d < dof - cdof; ++d) { 595cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 596cb1e1211SMatthew G Knepley 5979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(leafSectionAdj, off + d, &ldof)); 5989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(rootSectionAdj, off + d, &rdof)); 599cb1e1211SMatthew G Knepley if (ldof > 0) { 600cb1e1211SMatthew G Knepley /* We do not own this point */ 601cb1e1211SMatthew G Knepley } else if (rdof > 0) { 602cb1e1211SMatthew G Knepley PetscInt aoff, roff; 603cb1e1211SMatthew G Knepley 6049566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, goff + d, &aoff)); 6059566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(rootSectionAdj, off + d, &roff)); 6069566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(&cols[aoff], &rootAdj[roff], rdof)); 607cb1e1211SMatthew G Knepley } else { 608cb1e1211SMatthew G Knepley found = PETSC_FALSE; 609cb1e1211SMatthew G Knepley } 610cb1e1211SMatthew G Knepley } 611cb1e1211SMatthew G Knepley if (found) continue; 6129566063dSJacob Faibussowitsch PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj)); 613a38eeca9SJames Wright PetscCall(AdjancencyContainsLeafRootPair(numMyRankPair, leavesMyRankPair, rootsMyRankPair, numAdj, tmpAdj, &num_exclude_leaves, exclude_leaves)); 6149566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(anchorSectionAdj, p, &anDof)); 6159566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(anchorSectionAdj, p, &anOff)); 616cb1e1211SMatthew G Knepley for (d = goff; d < goff + dof - cdof; ++d) { 617cb1e1211SMatthew G Knepley PetscInt adof, aoff, i = 0; 618cb1e1211SMatthew G Knepley 6199566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, d, &adof)); 6209566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, d, &aoff)); 621cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 622a38eeca9SJames Wright const PetscInt padj = tmpAdj[q], *ncind; 623a38eeca9SJames Wright PetscInt ndof, ncdof, ngoff, nd, count; 624cb1e1211SMatthew G Knepley 625cb1e1211SMatthew G Knepley /* Adjacent points may not be in the section chart */ 626cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 6279566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(section, padj, &ndof)); 6289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintDof(section, padj, &ncdof)); 6299566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstraintIndices(section, padj, &ncind)); 6309566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionGlobal, padj, &ngoff)); 631a38eeca9SJames Wright // If leaf-root pair are both on this rank, only count root 632a38eeca9SJames Wright PetscCall(PetscFindInt(padj, num_exclude_leaves, exclude_leaves, &count)); 633a38eeca9SJames Wright if (count >= 0) continue; 634ad540459SPierre Jolivet for (nd = 0; nd < ndof - ncdof; ++nd, ++i) cols[aoff + i] = ngoff < 0 ? -(ngoff + 1) + nd : ngoff + nd; 635cb1e1211SMatthew G Knepley } 636ad540459SPierre Jolivet for (q = 0; q < anDof; q++, i++) cols[aoff + i] = anchorAdj[anOff + q]; 63763a3b9bcSJacob Faibussowitsch PetscCheck(i == adof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid number of entries %" PetscInt_FMT " != %" PetscInt_FMT " for dof %" PetscInt_FMT " (point %" PetscInt_FMT ")", i, adof, d, p); 638cb1e1211SMatthew G Knepley } 639cb1e1211SMatthew G Knepley } 6409566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&anchorSectionAdj)); 6419566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&leafSectionAdj)); 6429566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&rootSectionAdj)); 643a38eeca9SJames Wright PetscCall(PetscFree(exclude_leaves)); 644a38eeca9SJames Wright PetscCall(PetscFree2(rootsMyRankPair, leavesMyRankPair)); 6459566063dSJacob Faibussowitsch PetscCall(PetscFree(anchorAdj)); 6469566063dSJacob Faibussowitsch PetscCall(PetscFree(rootAdj)); 6479566063dSJacob Faibussowitsch PetscCall(PetscFree(tmpAdj)); 648cb1e1211SMatthew G Knepley /* Debugging */ 649cb1e1211SMatthew G Knepley if (debug) { 6509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Column indices\n")); 6516a1a2f7bSJames Wright PetscCall(PetscSectionArrayView(sectionAdj, cols, PETSC_INT, NULL)); 652cb1e1211SMatthew G Knepley } 653a9fb6443SMatthew G. Knepley 654a9fb6443SMatthew G. Knepley *sA = sectionAdj; 655a9fb6443SMatthew G. Knepley *colIdx = cols; 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 657a9fb6443SMatthew G. Knepley } 658a9fb6443SMatthew G. Knepley 659d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexUpdateAllocation_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[]) 660d71ae5a4SJacob Faibussowitsch { 661e101f074SMatthew G. Knepley PetscSection section; 662a9fb6443SMatthew G. Knepley PetscInt rStart, rEnd, r, pStart, pEnd, p; 663a9fb6443SMatthew G. Knepley 664a9fb6443SMatthew G. Knepley PetscFunctionBegin; 665a9fb6443SMatthew G. Knepley /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */ 6669566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd)); 6671dca8a05SBarry Smith PetscCheck((rStart % bs) == 0 && (rEnd % bs) == 0, PetscObjectComm((PetscObject)rLayout), PETSC_ERR_ARG_WRONG, "Invalid layout [%" PetscInt_FMT ", %" PetscInt_FMT ") for matrix, must be divisible by block size %" PetscInt_FMT, rStart, rEnd, bs); 668a9fb6443SMatthew G. Knepley if (f >= 0 && bs == 1) { 6699566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 6709566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 671a9fb6443SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 672294b7009SMatthew G. Knepley PetscInt rS, rE; 673a9fb6443SMatthew G. Knepley 6749566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE)); 675a9fb6443SMatthew G. Knepley for (r = rS; r < rE; ++r) { 676a9fb6443SMatthew G. Knepley PetscInt numCols, cStart, c; 677a9fb6443SMatthew G. Knepley 6789566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols)); 6799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart)); 680a9fb6443SMatthew G. Knepley for (c = cStart; c < cStart + numCols; ++c) { 681a9fb6443SMatthew G. Knepley if ((cols[c] >= rStart) && (cols[c] < rEnd)) { 682a9fb6443SMatthew G. Knepley ++dnz[r - rStart]; 683a9fb6443SMatthew G. Knepley if (cols[c] >= r) ++dnzu[r - rStart]; 684a9fb6443SMatthew G. Knepley } else { 685a9fb6443SMatthew G. Knepley ++onz[r - rStart]; 686a9fb6443SMatthew G. Knepley if (cols[c] >= r) ++onzu[r - rStart]; 687a9fb6443SMatthew G. Knepley } 688a9fb6443SMatthew G. Knepley } 689a9fb6443SMatthew G. Knepley } 690a9fb6443SMatthew G. Knepley } 691a9fb6443SMatthew G. Knepley } else { 692cb1e1211SMatthew G Knepley /* Only loop over blocks of rows */ 693cb1e1211SMatthew G Knepley for (r = rStart / bs; r < rEnd / bs; ++r) { 694cb1e1211SMatthew G Knepley const PetscInt row = r * bs; 695cb1e1211SMatthew G Knepley PetscInt numCols, cStart, c; 696cb1e1211SMatthew G Knepley 6979566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, row, &numCols)); 6989566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, row, &cStart)); 699cb1e1211SMatthew G Knepley for (c = cStart; c < cStart + numCols; ++c) { 700e7bcfa36SMatthew G. Knepley if ((cols[c] >= rStart) && (cols[c] < rEnd)) { 701e7bcfa36SMatthew G. Knepley ++dnz[r - rStart / bs]; 702e7bcfa36SMatthew G. Knepley if (cols[c] >= row) ++dnzu[r - rStart / bs]; 703cb1e1211SMatthew G Knepley } else { 704e7bcfa36SMatthew G. Knepley ++onz[r - rStart / bs]; 705e7bcfa36SMatthew G. Knepley if (cols[c] >= row) ++onzu[r - rStart / bs]; 706cb1e1211SMatthew G Knepley } 707cb1e1211SMatthew G Knepley } 708cb1e1211SMatthew G Knepley } 709a9fb6443SMatthew G. Knepley for (r = 0; r < (rEnd - rStart) / bs; ++r) { 710cb1e1211SMatthew G Knepley dnz[r] /= bs; 711cb1e1211SMatthew G Knepley onz[r] /= bs; 712cb1e1211SMatthew G Knepley dnzu[r] /= bs; 713cb1e1211SMatthew G Knepley onzu[r] /= bs; 714cb1e1211SMatthew G Knepley } 715cb1e1211SMatthew G Knepley } 7163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 717a9fb6443SMatthew G. Knepley } 718a9fb6443SMatthew G. Knepley 719d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A) 720d71ae5a4SJacob Faibussowitsch { 721e101f074SMatthew G. Knepley PetscSection section; 722a9fb6443SMatthew G. Knepley PetscScalar *values; 723a9fb6443SMatthew G. Knepley PetscInt rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0; 724a9fb6443SMatthew G. Knepley 725a9fb6443SMatthew G. Knepley PetscFunctionBegin; 7269566063dSJacob Faibussowitsch PetscCall(PetscLayoutGetRange(rLayout, &rStart, &rEnd)); 727a9fb6443SMatthew G. Knepley for (r = rStart; r < rEnd; ++r) { 7289566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &len)); 729a9fb6443SMatthew G. Knepley maxRowLen = PetscMax(maxRowLen, len); 730a9fb6443SMatthew G. Knepley } 7319566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(maxRowLen, &values)); 732a9fb6443SMatthew G. Knepley if (f >= 0 && bs == 1) { 7339566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 7349566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 735a9fb6443SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 736294b7009SMatthew G. Knepley PetscInt rS, rE; 737a9fb6443SMatthew G. Knepley 7389566063dSJacob Faibussowitsch PetscCall(DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE)); 739a9fb6443SMatthew G. Knepley for (r = rS; r < rE; ++r) { 7407e01ee02SMatthew G. Knepley PetscInt numCols, cStart; 741a9fb6443SMatthew G. Knepley 7429566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols)); 7439566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart)); 7449566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES)); 745a9fb6443SMatthew G. Knepley } 746a9fb6443SMatthew G. Knepley } 747a9fb6443SMatthew G. Knepley } else { 748a9fb6443SMatthew G. Knepley for (r = rStart; r < rEnd; ++r) { 749a9fb6443SMatthew G. Knepley PetscInt numCols, cStart; 750a9fb6443SMatthew G. Knepley 7519566063dSJacob Faibussowitsch PetscCall(PetscSectionGetDof(sectionAdj, r, &numCols)); 7529566063dSJacob Faibussowitsch PetscCall(PetscSectionGetOffset(sectionAdj, r, &cStart)); 7539566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES)); 754a9fb6443SMatthew G. Knepley } 755a9fb6443SMatthew G. Knepley } 7569566063dSJacob Faibussowitsch PetscCall(PetscFree(values)); 7573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 758a9fb6443SMatthew G. Knepley } 759a9fb6443SMatthew G. Knepley 760cc4c1da9SBarry Smith /*@ 761a1cb98faSBarry Smith DMPlexPreallocateOperator - Calculate the matrix nonzero pattern based upon the information in the `DM`, 762a1cb98faSBarry Smith the `PetscDS` it contains, and the default `PetscSection`. 76325e402d2SMatthew G. Knepley 76425e402d2SMatthew G. Knepley Collective 76525e402d2SMatthew G. Knepley 7664165533cSJose E. Roman Input Parameters: 767a1cb98faSBarry Smith + dm - The `DMPLEX` 76825e402d2SMatthew G. Knepley . bs - The matrix blocksize 76925e402d2SMatthew G. Knepley . dnz - An array to hold the number of nonzeros in the diagonal block 77025e402d2SMatthew G. Knepley . onz - An array to hold the number of nonzeros in the off-diagonal block 77125e402d2SMatthew G. Knepley . dnzu - An array to hold the number of nonzeros in the upper triangle of the diagonal block 77225e402d2SMatthew G. Knepley . onzu - An array to hold the number of nonzeros in the upper triangle of the off-diagonal block 773a1cb98faSBarry Smith - fillMatrix - If `PETSC_TRUE`, fill the matrix with zeros 77425e402d2SMatthew G. Knepley 7754165533cSJose E. Roman Output Parameter: 77625e402d2SMatthew G. Knepley . A - The preallocated matrix 77725e402d2SMatthew G. Knepley 77825e402d2SMatthew G. Knepley Level: advanced 77925e402d2SMatthew G. Knepley 7801cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `DMPLEX`, `DMCreateMatrix()` 78125e402d2SMatthew G. Knepley @*/ 782d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix) 783d71ae5a4SJacob Faibussowitsch { 784a9fb6443SMatthew G. Knepley MPI_Comm comm; 785a9fb6443SMatthew G. Knepley MatType mtype; 786a9fb6443SMatthew G. Knepley PetscSF sf, sfDof; 787e101f074SMatthew G. Knepley PetscSection section; 788a9fb6443SMatthew G. Knepley PetscInt *remoteOffsets; 789a9fb6443SMatthew G. Knepley PetscSection sectionAdj[4] = {NULL, NULL, NULL, NULL}; 790a9fb6443SMatthew G. Knepley PetscInt *cols[4] = {NULL, NULL, NULL, NULL}; 791a9fb6443SMatthew G. Knepley PetscBool useCone, useClosure; 7927e01ee02SMatthew G. Knepley PetscInt Nf, f, idx, locRows; 793a9fb6443SMatthew G. Knepley PetscLayout rLayout; 794a9fb6443SMatthew G. Knepley PetscBool isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE; 795a9fb6443SMatthew G. Knepley 796a9fb6443SMatthew G. Knepley PetscFunctionBegin; 797a9fb6443SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 798064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(A, MAT_CLASSID, 7); 7994f572ea9SToby Isaac if (dnz) PetscAssertPointer(dnz, 3); 8004f572ea9SToby Isaac if (onz) PetscAssertPointer(onz, 4); 8014f572ea9SToby Isaac if (dnzu) PetscAssertPointer(dnzu, 5); 8024f572ea9SToby Isaac if (onzu) PetscAssertPointer(onzu, 6); 803a38eeca9SJames Wright PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf)); 8049566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 8059566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-dm_view_preallocation", &debug, NULL)); 8069566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)dm, &comm)); 8079566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(DMPLEX_Preallocate, dm, 0, 0, 0)); 808a9fb6443SMatthew G. Knepley /* Create dof SF based on point SF */ 809a9fb6443SMatthew G. Knepley if (debug) { 810e101f074SMatthew G. Knepley PetscSection section, sectionGlobal; 811a9fb6443SMatthew G. Knepley PetscSF sf; 812a9fb6443SMatthew G. Knepley 813a38eeca9SJames Wright PetscCall(DMGetIsoperiodicPointSF_Internal(dm, &sf)); 8149566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm, §ion)); 8159566063dSJacob Faibussowitsch PetscCall(DMGetGlobalSection(dm, §ionGlobal)); 8169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Input Section for Preallocation:\n")); 8179566063dSJacob Faibussowitsch PetscCall(PetscSectionView(section, NULL)); 8189566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Input Global Section for Preallocation:\n")); 8199566063dSJacob Faibussowitsch PetscCall(PetscSectionView(sectionGlobal, NULL)); 8209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Input SF for Preallocation:\n")); 8219566063dSJacob Faibussowitsch PetscCall(PetscSFView(sf, NULL)); 822a9fb6443SMatthew G. Knepley } 8239566063dSJacob Faibussowitsch PetscCall(PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets)); 8249566063dSJacob Faibussowitsch PetscCall(PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof)); 8259566063dSJacob Faibussowitsch PetscCall(PetscFree(remoteOffsets)); 8266a1a2f7bSJames Wright if (debug) { 8279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "Dof SF for Preallocation:\n")); 8289566063dSJacob Faibussowitsch PetscCall(PetscSFView(sfDof, NULL)); 829a9fb6443SMatthew G. Knepley } 830a9fb6443SMatthew G. Knepley /* Create allocation vectors from adjacency graph */ 8319566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(A, &locRows, NULL)); 8329566063dSJacob Faibussowitsch PetscCall(PetscLayoutCreate(comm, &rLayout)); 8339566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetLocalSize(rLayout, locRows)); 8349566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetBlockSize(rLayout, 1)); 8359566063dSJacob Faibussowitsch PetscCall(PetscLayoutSetUp(rLayout)); 836a9fb6443SMatthew G. Knepley /* There are 4 types of adjacency */ 8379566063dSJacob Faibussowitsch PetscCall(PetscSectionGetNumFields(section, &Nf)); 838a9fb6443SMatthew G. Knepley if (Nf < 1 || bs > 1) { 8399566063dSJacob Faibussowitsch PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure)); 840a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 8419566063dSJacob Faibussowitsch PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, §ionAdj[idx], &cols[idx])); 8429566063dSJacob Faibussowitsch PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu)); 843a9fb6443SMatthew G. Knepley } else { 844a9fb6443SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 8459566063dSJacob Faibussowitsch PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure)); 846a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 8479566063dSJacob Faibussowitsch if (!sectionAdj[idx]) PetscCall(DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, §ionAdj[idx], &cols[idx])); 8489566063dSJacob Faibussowitsch PetscCall(DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu)); 849a9fb6443SMatthew G. Knepley } 850a9fb6443SMatthew G. Knepley } 8519566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&sfDof)); 852cb1e1211SMatthew G Knepley /* Set matrix pattern */ 8539566063dSJacob Faibussowitsch PetscCall(MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu)); 8549566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_TRUE)); 85589545effSMatthew G. Knepley /* Check for symmetric storage */ 8569566063dSJacob Faibussowitsch PetscCall(MatGetType(A, &mtype)); 8579566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype, MATSBAIJ, &isSymBlock)); 8589566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock)); 8599566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock)); 8609566063dSJacob Faibussowitsch if (isSymBlock || isSymSeqBlock || isSymMPIBlock) PetscCall(MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE)); 861cb1e1211SMatthew G Knepley /* Fill matrix with zeros */ 862cb1e1211SMatthew G Knepley if (fillMatrix) { 863a9fb6443SMatthew G. Knepley if (Nf < 1 || bs > 1) { 8649566063dSJacob Faibussowitsch PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure)); 865a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 8669566063dSJacob Faibussowitsch PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A)); 867a9fb6443SMatthew G. Knepley } else { 868a9fb6443SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 8699566063dSJacob Faibussowitsch PetscCall(DMGetAdjacency(dm, f, &useCone, &useClosure)); 870a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 8719566063dSJacob Faibussowitsch PetscCall(DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A)); 872cb1e1211SMatthew G Knepley } 873cb1e1211SMatthew G Knepley } 8749566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 8759566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 876cb1e1211SMatthew G Knepley } 8779566063dSJacob Faibussowitsch PetscCall(PetscLayoutDestroy(&rLayout)); 8789371c9d4SSatish Balay for (idx = 0; idx < 4; ++idx) { 8799371c9d4SSatish Balay PetscCall(PetscSectionDestroy(§ionAdj[idx])); 8809371c9d4SSatish Balay PetscCall(PetscFree(cols[idx])); 8819371c9d4SSatish Balay } 8829566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(DMPLEX_Preallocate, dm, 0, 0, 0)); 8833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 884cb1e1211SMatthew G Knepley } 885cb1e1211SMatthew G Knepley 886cb1e1211SMatthew G Knepley #if 0 887cb1e1211SMatthew 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) 888cb1e1211SMatthew G Knepley { 889cb1e1211SMatthew G Knepley PetscInt *tmpClosure,*tmpAdj,*visits; 890cb1e1211SMatthew G Knepley PetscInt c,cStart,cEnd,pStart,pEnd; 891cb1e1211SMatthew G Knepley 892cb1e1211SMatthew G Knepley PetscFunctionBegin; 8939566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 8949566063dSJacob Faibussowitsch PetscCall(DMPlexGetDepth(dm, &depth)); 8959566063dSJacob Faibussowitsch PetscCall(DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize)); 896cb1e1211SMatthew G Knepley 897e7b80042SMatthew G. Knepley maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1)); 898cb1e1211SMatthew G Knepley 8999566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 900cb1e1211SMatthew G Knepley npoints = pEnd - pStart; 901cb1e1211SMatthew G Knepley 9029566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits)); 9039566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(lvisits,pEnd-pStart)); 9049566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(visits,pEnd-pStart)); 9059566063dSJacob Faibussowitsch PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 906cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 907cb1e1211SMatthew G Knepley PetscInt *support = tmpClosure; 9089566063dSJacob Faibussowitsch PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support)); 909cb1e1211SMatthew G Knepley for (p=0; p<supportSize; p++) lvisits[support[p]]++; 910cb1e1211SMatthew G Knepley } 9119566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM)); 9129566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd (sf,MPIU_INT,lvisits,visits,MPI_SUM)); 9139566063dSJacob Faibussowitsch PetscCall(PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE)); 9149566063dSJacob Faibussowitsch PetscCall(PetscSFBcastEnd (sf,MPIU_INT,visits,lvisits)); 915cb1e1211SMatthew G Knepley 9169566063dSJacob Faibussowitsch PetscCall(PetscSFGetRootRanks()); 917cb1e1211SMatthew G Knepley 9189566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner)); 919cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 9209566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(cellmat,maxClosureSize*maxClosureSize)); 921cb1e1211SMatthew G Knepley /* 922cb1e1211SMatthew G Knepley Depth-first walk of transitive closure. 923cb1e1211SMatthew 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. 924cb1e1211SMatthew G Knepley This contribution is added to dnz if owning ranks of p and q match, to onz otherwise. 925cb1e1211SMatthew G Knepley */ 926cb1e1211SMatthew G Knepley } 927cb1e1211SMatthew G Knepley 9289566063dSJacob Faibussowitsch PetscCall(PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM)); 9299566063dSJacob Faibussowitsch PetscCall(PetscSFReduceEnd (sf,MPIU_INT,lonz,onz,MPI_SUM)); 9303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 931cb1e1211SMatthew G Knepley } 932cb1e1211SMatthew G Knepley #endif 933