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() */ 7e101f074SMatthew G. Knepley static PetscErrorCode DMPlexComputeAnchorAdjacencies(DM dm, PetscBool useCone, PetscBool useClosure, PetscSection *anchorSectionAdj, PetscInt *anchorAdj[]) 876185916SToby Isaac { 976185916SToby Isaac PetscInt pStart, pEnd; 10e101f074SMatthew G. Knepley PetscSection section, sectionGlobal, adjSec, aSec; 1176185916SToby Isaac IS aIS; 1276185916SToby Isaac PetscErrorCode ierr; 1376185916SToby Isaac 1476185916SToby Isaac PetscFunctionBegin; 1592fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 16e87a4003SBarry Smith ierr = DMGetGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 1776185916SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject) section), &adjSec);CHKERRQ(ierr); 1876185916SToby Isaac ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 1976185916SToby Isaac ierr = PetscSectionSetChart(adjSec,pStart,pEnd);CHKERRQ(ierr); 2076185916SToby Isaac 21a17985deSToby Isaac ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr); 2276185916SToby Isaac if (aSec) { 2376185916SToby Isaac const PetscInt *anchors; 2476185916SToby Isaac PetscInt p, q, a, aSize, *offsets, aStart, aEnd, *inverse, iSize, *adj, adjSize; 2576185916SToby Isaac PetscInt *tmpAdjP = NULL, *tmpAdjQ = NULL; 2676185916SToby Isaac PetscSection inverseSec; 2776185916SToby Isaac 2876185916SToby Isaac /* invert the constraint-to-anchor map */ 2976185916SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)aSec),&inverseSec);CHKERRQ(ierr); 3076185916SToby Isaac ierr = PetscSectionSetChart(inverseSec,pStart,pEnd);CHKERRQ(ierr); 3176185916SToby Isaac ierr = ISGetLocalSize(aIS, &aSize);CHKERRQ(ierr); 3276185916SToby Isaac ierr = ISGetIndices(aIS, &anchors);CHKERRQ(ierr); 3376185916SToby Isaac 3476185916SToby Isaac for (p = 0; p < aSize; p++) { 3576185916SToby Isaac PetscInt a = anchors[p]; 3676185916SToby Isaac 3776185916SToby Isaac ierr = PetscSectionAddDof(inverseSec,a,1);CHKERRQ(ierr); 3876185916SToby Isaac } 3976185916SToby Isaac ierr = PetscSectionSetUp(inverseSec);CHKERRQ(ierr); 4076185916SToby Isaac ierr = PetscSectionGetStorageSize(inverseSec,&iSize);CHKERRQ(ierr); 4176185916SToby Isaac ierr = PetscMalloc1(iSize,&inverse);CHKERRQ(ierr); 4276185916SToby Isaac ierr = PetscCalloc1(pEnd-pStart,&offsets);CHKERRQ(ierr); 4376185916SToby Isaac ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr); 4476185916SToby Isaac for (p = aStart; p < aEnd; p++) { 4576185916SToby Isaac PetscInt dof, off; 4676185916SToby Isaac 4776185916SToby Isaac ierr = PetscSectionGetDof(aSec, p, &dof);CHKERRQ(ierr); 4876185916SToby Isaac ierr = PetscSectionGetOffset(aSec, p, &off);CHKERRQ(ierr); 4976185916SToby Isaac 5076185916SToby Isaac for (q = 0; q < dof; q++) { 5176185916SToby Isaac PetscInt iOff; 5276185916SToby Isaac 5376185916SToby Isaac a = anchors[off + q]; 5476185916SToby Isaac ierr = PetscSectionGetOffset(inverseSec, a, &iOff);CHKERRQ(ierr); 5576185916SToby Isaac inverse[iOff + offsets[a-pStart]++] = p; 5676185916SToby Isaac } 5776185916SToby Isaac } 5876185916SToby Isaac ierr = ISRestoreIndices(aIS, &anchors);CHKERRQ(ierr); 5976185916SToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 6076185916SToby Isaac 6176185916SToby Isaac /* construct anchorAdj and adjSec 6276185916SToby Isaac * 6376185916SToby Isaac * loop over anchors: 6476185916SToby Isaac * construct anchor adjacency 6576185916SToby Isaac * loop over constrained: 6676185916SToby Isaac * construct constrained adjacency 6776185916SToby Isaac * if not in anchor adjacency, add to dofs 6876185916SToby Isaac * setup adjSec, allocate anchorAdj 6976185916SToby Isaac * loop over anchors: 7076185916SToby Isaac * construct anchor adjacency 7176185916SToby Isaac * loop over constrained: 7276185916SToby Isaac * construct constrained adjacency 7376185916SToby Isaac * if not in anchor adjacency 7476185916SToby Isaac * if not already in list, put in list 7576185916SToby Isaac * sort, unique, reduce dof count 7676185916SToby Isaac * optional: compactify 7776185916SToby Isaac */ 7876185916SToby Isaac for (p = pStart; p < pEnd; p++) { 7976185916SToby Isaac PetscInt iDof, iOff, i, r, s, numAdjP = PETSC_DETERMINE; 8076185916SToby Isaac 8176185916SToby Isaac ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr); 8276185916SToby Isaac if (!iDof) continue; 8376185916SToby Isaac ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr); 84a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP);CHKERRQ(ierr); 8576185916SToby Isaac for (i = 0; i < iDof; i++) { 8676185916SToby Isaac PetscInt iNew = 0, qAdj, qAdjDof, qAdjCDof, numAdjQ = PETSC_DETERMINE; 8776185916SToby Isaac 8876185916SToby Isaac q = inverse[iOff + i]; 89a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm,q,useCone,useClosure,PETSC_TRUE,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr); 9076185916SToby Isaac for (r = 0; r < numAdjQ; r++) { 9176185916SToby Isaac qAdj = tmpAdjQ[r]; 9276185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue; 9376185916SToby Isaac for (s = 0; s < numAdjP; s++) { 9476185916SToby Isaac if (qAdj == tmpAdjP[s]) break; 9576185916SToby Isaac } 9676185916SToby Isaac if (s < numAdjP) continue; 9776185916SToby Isaac ierr = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr); 9876185916SToby Isaac ierr = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr); 9976185916SToby Isaac iNew += qAdjDof - qAdjCDof; 10076185916SToby Isaac } 10176185916SToby Isaac ierr = PetscSectionAddDof(adjSec,p,iNew);CHKERRQ(ierr); 10276185916SToby Isaac } 10376185916SToby Isaac } 10476185916SToby Isaac 10576185916SToby Isaac ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr); 10676185916SToby Isaac ierr = PetscSectionGetStorageSize(adjSec,&adjSize);CHKERRQ(ierr); 10776185916SToby Isaac ierr = PetscMalloc1(adjSize,&adj);CHKERRQ(ierr); 10876185916SToby Isaac 10976185916SToby Isaac for (p = pStart; p < pEnd; p++) { 11076185916SToby Isaac PetscInt iDof, iOff, i, r, s, aOff, aOffOrig, aDof, numAdjP = PETSC_DETERMINE; 11176185916SToby Isaac 11276185916SToby Isaac ierr = PetscSectionGetDof(inverseSec,p,&iDof);CHKERRQ(ierr); 11376185916SToby Isaac if (!iDof) continue; 11476185916SToby Isaac ierr = PetscSectionGetOffset(inverseSec,p,&iOff);CHKERRQ(ierr); 115a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm,p,useCone,useClosure,PETSC_TRUE,&numAdjP,&tmpAdjP);CHKERRQ(ierr); 11676185916SToby Isaac ierr = PetscSectionGetDof(adjSec,p,&aDof);CHKERRQ(ierr); 11776185916SToby Isaac ierr = PetscSectionGetOffset(adjSec,p,&aOff);CHKERRQ(ierr); 11876185916SToby Isaac aOffOrig = aOff; 11976185916SToby Isaac for (i = 0; i < iDof; i++) { 12076185916SToby Isaac PetscInt qAdj, qAdjDof, qAdjCDof, qAdjOff, nd, numAdjQ = PETSC_DETERMINE; 12176185916SToby Isaac 12276185916SToby Isaac q = inverse[iOff + i]; 123a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm,q,useCone,useClosure,PETSC_TRUE,&numAdjQ,&tmpAdjQ);CHKERRQ(ierr); 12476185916SToby Isaac for (r = 0; r < numAdjQ; r++) { 12576185916SToby Isaac qAdj = tmpAdjQ[r]; 12676185916SToby Isaac if ((qAdj < pStart) || (qAdj >= pEnd)) continue; 12776185916SToby Isaac for (s = 0; s < numAdjP; s++) { 12876185916SToby Isaac if (qAdj == tmpAdjP[s]) break; 12976185916SToby Isaac } 13076185916SToby Isaac if (s < numAdjP) continue; 13176185916SToby Isaac ierr = PetscSectionGetDof(section,qAdj,&qAdjDof);CHKERRQ(ierr); 13276185916SToby Isaac ierr = PetscSectionGetConstraintDof(section,qAdj,&qAdjCDof);CHKERRQ(ierr); 13376185916SToby Isaac ierr = PetscSectionGetOffset(sectionGlobal,qAdj,&qAdjOff);CHKERRQ(ierr); 13476185916SToby Isaac for (nd = 0; nd < qAdjDof-qAdjCDof; ++nd) { 13576185916SToby Isaac adj[aOff++] = (qAdjOff < 0 ? -(qAdjOff+1) : qAdjOff) + nd; 13676185916SToby Isaac } 13776185916SToby Isaac } 13876185916SToby Isaac } 13976185916SToby Isaac ierr = PetscSortRemoveDupsInt(&aDof,&adj[aOffOrig]);CHKERRQ(ierr); 14076185916SToby Isaac ierr = PetscSectionSetDof(adjSec,p,aDof);CHKERRQ(ierr); 14176185916SToby Isaac } 14276185916SToby Isaac *anchorAdj = adj; 14376185916SToby Isaac 14476185916SToby Isaac /* clean up */ 14576185916SToby Isaac ierr = PetscSectionDestroy(&inverseSec);CHKERRQ(ierr); 14676185916SToby Isaac ierr = PetscFree(inverse);CHKERRQ(ierr); 14776185916SToby Isaac ierr = PetscFree(tmpAdjP);CHKERRQ(ierr); 14876185916SToby Isaac ierr = PetscFree(tmpAdjQ);CHKERRQ(ierr); 14976185916SToby Isaac } 15076185916SToby Isaac else { 15176185916SToby Isaac *anchorAdj = NULL; 15276185916SToby Isaac ierr = PetscSectionSetUp(adjSec);CHKERRQ(ierr); 15376185916SToby Isaac } 15476185916SToby Isaac *anchorSectionAdj = adjSec; 15576185916SToby Isaac PetscFunctionReturn(0); 15676185916SToby Isaac } 15776185916SToby Isaac 15825e402d2SMatthew G. Knepley static PetscErrorCode DMPlexCreateAdjacencySection_Static(DM dm, PetscInt bs, PetscSF sfDof, PetscBool useCone, PetscBool useClosure, PetscBool useAnchors, PetscSection *sA, PetscInt **colIdx) 159cb1e1211SMatthew G Knepley { 160cb1e1211SMatthew G Knepley MPI_Comm comm; 161a9fb6443SMatthew G. Knepley PetscMPIInt size; 162a9fb6443SMatthew G. Knepley PetscBool doCommLocal, doComm, debug = PETSC_FALSE; 163a9fb6443SMatthew G. Knepley PetscSF sf, sfAdj; 164e101f074SMatthew G. Knepley PetscSection section, sectionGlobal, leafSectionAdj, rootSectionAdj, sectionAdj, anchorSectionAdj; 165a9fb6443SMatthew G. Knepley PetscInt nroots, nleaves, l, p, r; 166cb1e1211SMatthew G Knepley const PetscInt *leaves; 167cb1e1211SMatthew G Knepley const PetscSFNode *remotes; 168cb1e1211SMatthew G Knepley PetscInt dim, pStart, pEnd, numDof, globalOffStart, globalOffEnd, numCols; 1698821704fSMatthew G. Knepley PetscInt *tmpAdj = NULL, *adj, *rootAdj, *anchorAdj = NULL, *cols, *remoteOffsets; 17070034214SMatthew G. Knepley PetscInt adjSize; 171cb1e1211SMatthew G Knepley PetscErrorCode ierr; 172cb1e1211SMatthew G Knepley 173cb1e1211SMatthew G Knepley PetscFunctionBegin; 174cb1e1211SMatthew G Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 175c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr); 176ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 177c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 178cb1e1211SMatthew G Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 17992fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 180e87a4003SBarry Smith ierr = DMGetGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 18147a6104aSMatthew G. Knepley ierr = PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 18247a6104aSMatthew G. Knepley doCommLocal = (size > 1) && (nroots >= 0) ? PETSC_TRUE : PETSC_FALSE; 183820f2d46SBarry Smith ierr = MPIU_Allreduce(&doCommLocal, &doComm, 1, MPIU_BOOL, MPI_LAND, comm);CHKERRMPI(ierr); 184cb1e1211SMatthew G Knepley /* Create section for dof adjacency (dof ==> # adj dof) */ 185cb1e1211SMatthew G Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 186cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(section, &numDof);CHKERRQ(ierr); 187cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, &leafSectionAdj);CHKERRQ(ierr); 188cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(leafSectionAdj, 0, numDof);CHKERRQ(ierr); 189cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, &rootSectionAdj);CHKERRQ(ierr); 190cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(rootSectionAdj, 0, numDof);CHKERRQ(ierr); 191cb1e1211SMatthew G Knepley /* Fill in the ghost dofs on the interface */ 192cb1e1211SMatthew G Knepley ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, &remotes);CHKERRQ(ierr); 193cb1e1211SMatthew G Knepley /* 194964bf7afSMatthew G. Knepley section - maps points to (# dofs, local dofs) 195964bf7afSMatthew G. Knepley sectionGlobal - maps points to (# dofs, global dofs) 196964bf7afSMatthew G. Knepley leafSectionAdj - maps unowned local dofs to # adj dofs 197964bf7afSMatthew G. Knepley rootSectionAdj - maps owned local dofs to # adj dofs 198964bf7afSMatthew G. Knepley adj - adj global dofs indexed by leafSectionAdj 199964bf7afSMatthew G. Knepley rootAdj - adj global dofs indexed by rootSectionAdj 200964bf7afSMatthew G. Knepley sf - describes shared points across procs 201964bf7afSMatthew G. Knepley sfDof - describes shared dofs across procs 202964bf7afSMatthew G. Knepley sfAdj - describes shared adjacent dofs across procs 203cb1e1211SMatthew G Knepley ** The bootstrapping process involves six rounds with similar structure of visiting neighbors of each point. 20476185916SToby Isaac (0). If there are point-to-point constraints, add the adjacencies of constrained points to anchors in anchorAdj 20576185916SToby Isaac (This is done in DMPlexComputeAnchorAdjacencies()) 206cb1e1211SMatthew G Knepley 1. Visit unowned points on interface, count adjacencies placing in leafSectionAdj 207cb1e1211SMatthew G Knepley Reduce those counts to rootSectionAdj (now redundantly counting some interface points) 208cb1e1211SMatthew G Knepley 2. Visit owned points on interface, count adjacencies placing in rootSectionAdj 209cb1e1211SMatthew G Knepley Create sfAdj connecting rootSectionAdj and leafSectionAdj 210cb1e1211SMatthew G Knepley 3. Visit unowned points on interface, write adjacencies to adj 211cb1e1211SMatthew G Knepley Gather adj to rootAdj (note that there is redundancy in rootAdj when multiple procs find the same adjacencies) 212cb1e1211SMatthew G Knepley 4. Visit owned points on interface, write adjacencies to rootAdj 213cb1e1211SMatthew G Knepley Remove redundancy in rootAdj 214cb1e1211SMatthew G Knepley ** The last two traversals use transitive closure 215cb1e1211SMatthew G Knepley 5. Visit all owned points in the subdomain, count dofs for each point (sectionAdj) 216cb1e1211SMatthew G Knepley Allocate memory addressed by sectionAdj (cols) 217cb1e1211SMatthew G Knepley 6. Visit all owned points in the subdomain, insert dof adjacencies into cols 218cb1e1211SMatthew G Knepley ** Knowing all the column adjacencies, check ownership and sum into dnz and onz 219cb1e1211SMatthew G Knepley */ 220e101f074SMatthew G. Knepley ierr = DMPlexComputeAnchorAdjacencies(dm, useCone, useClosure, &anchorSectionAdj, &anchorAdj);CHKERRQ(ierr); 221cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 22276185916SToby Isaac PetscInt dof, off, d, q, anDof; 22370034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 224cb1e1211SMatthew G Knepley 225fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 226cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 227cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 228a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr); 229cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 230cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 231cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 232cb1e1211SMatthew G Knepley 233cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 234cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 235cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 236cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 237cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(leafSectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 238cb1e1211SMatthew G Knepley } 239cb1e1211SMatthew G Knepley } 24076185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 24176185916SToby Isaac if (anDof) { 24276185916SToby Isaac for (d = off; d < off+dof; ++d) { 24376185916SToby Isaac ierr = PetscSectionAddDof(leafSectionAdj, d, anDof);CHKERRQ(ierr); 24476185916SToby Isaac } 24576185916SToby Isaac } 246cb1e1211SMatthew G Knepley } 247cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr); 248cb1e1211SMatthew G Knepley if (debug) { 249cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency Section for Preallocation on Leaves:\n");CHKERRQ(ierr); 250a8c5bd36SStefano Zampini ierr = PetscSectionView(leafSectionAdj, NULL);CHKERRQ(ierr); 251cb1e1211SMatthew G Knepley } 252cb1e1211SMatthew G Knepley /* Get maximum remote adjacency sizes for owned dofs on interface (roots) */ 25347a6104aSMatthew G. Knepley if (doComm) { 254cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr); 255cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd(sfDof, MPIU_INT, leafSectionAdj->atlasDof, rootSectionAdj->atlasDof, MPI_SUM);CHKERRQ(ierr); 256cb1e1211SMatthew G Knepley } 257cb1e1211SMatthew G Knepley if (debug) { 258cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots:\n");CHKERRQ(ierr); 259a8c5bd36SStefano Zampini ierr = PetscSectionView(rootSectionAdj, NULL);CHKERRQ(ierr); 260cb1e1211SMatthew G Knepley } 261cb1e1211SMatthew G Knepley /* Add in local adjacency sizes for owned dofs on interface (roots) */ 262cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 26376185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof; 264cb1e1211SMatthew G Knepley 265cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 266cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 267cb1e1211SMatthew G Knepley if (!dof) continue; 268cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 269cb1e1211SMatthew G Knepley if (adof <= 0) continue; 270a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr); 271cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 272cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 273cb1e1211SMatthew G Knepley PetscInt ndof, ncdof; 274cb1e1211SMatthew G Knepley 275cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 276cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 277cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 278cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 279cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(rootSectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 280cb1e1211SMatthew G Knepley } 281cb1e1211SMatthew G Knepley } 28276185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 28376185916SToby Isaac if (anDof) { 28476185916SToby Isaac for (d = off; d < off+dof; ++d) { 28576185916SToby Isaac ierr = PetscSectionAddDof(rootSectionAdj, d, anDof);CHKERRQ(ierr); 28676185916SToby Isaac } 28776185916SToby Isaac } 288cb1e1211SMatthew G Knepley } 289cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr); 290cb1e1211SMatthew G Knepley if (debug) { 291cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after local additions:\n");CHKERRQ(ierr); 292a8c5bd36SStefano Zampini ierr = PetscSectionView(rootSectionAdj, NULL);CHKERRQ(ierr); 293cb1e1211SMatthew G Knepley } 294cb1e1211SMatthew G Knepley /* Create adj SF based on dof SF */ 295cb1e1211SMatthew G Knepley ierr = PetscSFCreateRemoteOffsets(sfDof, rootSectionAdj, leafSectionAdj, &remoteOffsets);CHKERRQ(ierr); 296cb1e1211SMatthew G Knepley ierr = PetscSFCreateSectionSF(sfDof, rootSectionAdj, remoteOffsets, leafSectionAdj, &sfAdj);CHKERRQ(ierr); 2978a713f3bSLawrence Mitchell ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 2986cded2eaSMatthew G. Knepley if (debug && size > 1) { 299cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency SF for Preallocation:\n");CHKERRQ(ierr); 300cb1e1211SMatthew G Knepley ierr = PetscSFView(sfAdj, NULL);CHKERRQ(ierr); 301cb1e1211SMatthew G Knepley } 302cb1e1211SMatthew G Knepley /* Create leaf adjacency */ 303cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(leafSectionAdj);CHKERRQ(ierr); 304cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(leafSectionAdj, &adjSize);CHKERRQ(ierr); 3051795a4d1SJed Brown ierr = PetscCalloc1(adjSize, &adj);CHKERRQ(ierr); 306cb1e1211SMatthew G Knepley for (l = 0; l < nleaves; ++l) { 30776185916SToby Isaac PetscInt dof, off, d, q, anDof, anOff; 30870034214SMatthew G. Knepley PetscInt p = leaves[l], numAdj = PETSC_DETERMINE; 309cb1e1211SMatthew G Knepley 310fb47e2feSMatthew G. Knepley if ((p < pStart) || (p >= pEnd)) continue; 311cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 312cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 313a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr); 31476185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 31576185916SToby Isaac ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr); 316cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 317cb1e1211SMatthew G Knepley PetscInt aoff, i = 0; 318cb1e1211SMatthew G Knepley 319cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(leafSectionAdj, d, &aoff);CHKERRQ(ierr); 320cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 321cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 322cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 323cb1e1211SMatthew G Knepley 324cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 325cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 326cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 327cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 328cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd) { 329cb1e1211SMatthew G Knepley adj[aoff+i] = (ngoff < 0 ? -(ngoff+1) : ngoff) + nd; 330cb1e1211SMatthew G Knepley ++i; 331cb1e1211SMatthew G Knepley } 332cb1e1211SMatthew G Knepley } 33376185916SToby Isaac for (q = 0; q < anDof; q++) { 33476185916SToby Isaac adj[aoff+i] = anchorAdj[anOff+q]; 33576185916SToby Isaac ++i; 33676185916SToby Isaac } 337cb1e1211SMatthew G Knepley } 338cb1e1211SMatthew G Knepley } 339cb1e1211SMatthew G Knepley /* Debugging */ 340cb1e1211SMatthew G Knepley if (debug) { 341cb1e1211SMatthew G Knepley IS tmp; 342cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Leaf adjacency indices\n");CHKERRQ(ierr); 343cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, adj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 344cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 3450a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 346cb1e1211SMatthew G Knepley } 347543482b8SMatthew G. Knepley /* Gather adjacent indices to root */ 348cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(rootSectionAdj, &adjSize);CHKERRQ(ierr); 349785e854fSJed Brown ierr = PetscMalloc1(adjSize, &rootAdj);CHKERRQ(ierr); 350cb1e1211SMatthew G Knepley for (r = 0; r < adjSize; ++r) rootAdj[r] = -1; 35147a6104aSMatthew G. Knepley if (doComm) { 352543482b8SMatthew G. Knepley const PetscInt *indegree; 353543482b8SMatthew G. Knepley PetscInt *remoteadj, radjsize = 0; 354543482b8SMatthew G. Knepley 355543482b8SMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(sfAdj, &indegree);CHKERRQ(ierr); 356543482b8SMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(sfAdj, &indegree);CHKERRQ(ierr); 357543482b8SMatthew G. Knepley for (p = 0; p < adjSize; ++p) radjsize += indegree[p]; 358543482b8SMatthew G. Knepley ierr = PetscMalloc1(radjsize, &remoteadj);CHKERRQ(ierr); 359543482b8SMatthew G. Knepley ierr = PetscSFGatherBegin(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr); 360543482b8SMatthew G. Knepley ierr = PetscSFGatherEnd(sfAdj, MPIU_INT, adj, remoteadj);CHKERRQ(ierr); 3616dba2905SMatthew G. Knepley for (p = 0, l = 0, r = 0; p < adjSize; ++p, l = PetscMax(p, l + indegree[p-1])) { 362543482b8SMatthew G. Knepley PetscInt s; 3636dba2905SMatthew G. Knepley for (s = 0; s < indegree[p]; ++s, ++r) rootAdj[l+s] = remoteadj[r]; 364543482b8SMatthew G. Knepley } 365543482b8SMatthew G. Knepley if (r != radjsize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", r, radjsize); 3666dba2905SMatthew G. Knepley if (l != adjSize) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Inconsistency in communication %d != %d", l, adjSize); 367543482b8SMatthew G. Knepley ierr = PetscFree(remoteadj);CHKERRQ(ierr); 368cb1e1211SMatthew G Knepley } 369cb1e1211SMatthew G Knepley ierr = PetscSFDestroy(&sfAdj);CHKERRQ(ierr); 370cb1e1211SMatthew G Knepley ierr = PetscFree(adj);CHKERRQ(ierr); 371cb1e1211SMatthew G Knepley /* Debugging */ 372cb1e1211SMatthew G Knepley if (debug) { 373cb1e1211SMatthew G Knepley IS tmp; 374cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices after gather\n");CHKERRQ(ierr); 375cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 376cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 3770a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 378cb1e1211SMatthew G Knepley } 379cb1e1211SMatthew G Knepley /* Add in local adjacency indices for owned dofs on interface (roots) */ 380cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 38176185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, adof, dof, off, d, q, anDof, anOff; 382cb1e1211SMatthew G Knepley 383cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 384cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 385cb1e1211SMatthew G Knepley if (!dof) continue; 386cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 387cb1e1211SMatthew G Knepley if (adof <= 0) continue; 388a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr); 38976185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 39076185916SToby Isaac ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr); 391cb1e1211SMatthew G Knepley for (d = off; d < off+dof; ++d) { 392cb1e1211SMatthew G Knepley PetscInt adof, aoff, i; 393cb1e1211SMatthew G Knepley 394cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr); 395cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr); 396cb1e1211SMatthew G Knepley i = adof-1; 39776185916SToby Isaac for (q = 0; q < anDof; q++) { 39876185916SToby Isaac rootAdj[aoff+i] = anchorAdj[anOff+q]; 39976185916SToby Isaac --i; 40076185916SToby Isaac } 401cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 402cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 403cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 404cb1e1211SMatthew G Knepley 405cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 406cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 407cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 408cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 409cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd) { 410cb1e1211SMatthew G Knepley rootAdj[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd; 411cb1e1211SMatthew G Knepley --i; 412cb1e1211SMatthew G Knepley } 413cb1e1211SMatthew G Knepley } 414cb1e1211SMatthew G Knepley } 415cb1e1211SMatthew G Knepley } 416cb1e1211SMatthew G Knepley /* Debugging */ 417cb1e1211SMatthew G Knepley if (debug) { 418cb1e1211SMatthew G Knepley IS tmp; 419cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices\n");CHKERRQ(ierr); 420cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 421cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 4220a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 423cb1e1211SMatthew G Knepley } 424cb1e1211SMatthew G Knepley /* Compress indices */ 425cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(rootSectionAdj);CHKERRQ(ierr); 426cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 427cb1e1211SMatthew G Knepley PetscInt dof, cdof, off, d; 428cb1e1211SMatthew G Knepley PetscInt adof, aoff; 429cb1e1211SMatthew G Knepley 430cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 431cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 432cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 433cb1e1211SMatthew G Knepley if (!dof) continue; 434cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off, &adof);CHKERRQ(ierr); 435cb1e1211SMatthew G Knepley if (adof <= 0) continue; 436cb1e1211SMatthew G Knepley for (d = off; d < off+dof-cdof; ++d) { 437cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, d, &adof);CHKERRQ(ierr); 438cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, d, &aoff);CHKERRQ(ierr); 439cb1e1211SMatthew G Knepley ierr = PetscSortRemoveDupsInt(&adof, &rootAdj[aoff]);CHKERRQ(ierr); 440cb1e1211SMatthew G Knepley ierr = PetscSectionSetDof(rootSectionAdj, d, adof);CHKERRQ(ierr); 441cb1e1211SMatthew G Knepley } 442cb1e1211SMatthew G Knepley } 443cb1e1211SMatthew G Knepley /* Debugging */ 444cb1e1211SMatthew G Knepley if (debug) { 445cb1e1211SMatthew G Knepley IS tmp; 446cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjancency Section for Preallocation on Roots after compression:\n");CHKERRQ(ierr); 447a8c5bd36SStefano Zampini ierr = PetscSectionView(rootSectionAdj, NULL);CHKERRQ(ierr); 448cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Root adjacency indices after compression\n");CHKERRQ(ierr); 449cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, adjSize, rootAdj, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 450cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 4510a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 452cb1e1211SMatthew G Knepley } 453cb1e1211SMatthew G Knepley /* Build adjacency section: Maps global indices to sets of adjacent global indices */ 454cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffsetRange(sectionGlobal, &globalOffStart, &globalOffEnd);CHKERRQ(ierr); 455cb1e1211SMatthew G Knepley ierr = PetscSectionCreate(comm, §ionAdj);CHKERRQ(ierr); 456cb1e1211SMatthew G Knepley ierr = PetscSectionSetChart(sectionAdj, globalOffStart, globalOffEnd);CHKERRQ(ierr); 457cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 45876185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof; 459cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 460cb1e1211SMatthew G Knepley 461cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 462cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 463cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 464cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 465cb1e1211SMatthew G Knepley for (d = 0; d < dof-cdof; ++d) { 466cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 467cb1e1211SMatthew G Knepley 468cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr); 469cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr); 470cb1e1211SMatthew G Knepley if (ldof > 0) { 471cb1e1211SMatthew G Knepley /* We do not own this point */ 472cb1e1211SMatthew G Knepley } else if (rdof > 0) { 473cb1e1211SMatthew G Knepley ierr = PetscSectionSetDof(sectionAdj, goff+d, rdof);CHKERRQ(ierr); 474cb1e1211SMatthew G Knepley } else { 475cb1e1211SMatthew G Knepley found = PETSC_FALSE; 476cb1e1211SMatthew G Knepley } 477cb1e1211SMatthew G Knepley } 478cb1e1211SMatthew G Knepley if (found) continue; 479cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 480cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 481a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr); 482cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 483cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 484cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, noff; 485cb1e1211SMatthew G Knepley 486cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 487cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 488cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 489cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, padj, &noff);CHKERRQ(ierr); 490cb1e1211SMatthew G Knepley for (d = goff; d < goff+dof-cdof; ++d) { 491cb1e1211SMatthew G Knepley ierr = PetscSectionAddDof(sectionAdj, d, ndof-ncdof);CHKERRQ(ierr); 492cb1e1211SMatthew G Knepley } 493cb1e1211SMatthew G Knepley } 49476185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 49576185916SToby Isaac if (anDof) { 49676185916SToby Isaac for (d = goff; d < goff+dof-cdof; ++d) { 49776185916SToby Isaac ierr = PetscSectionAddDof(sectionAdj, d, anDof);CHKERRQ(ierr); 49876185916SToby Isaac } 49976185916SToby Isaac } 500cb1e1211SMatthew G Knepley } 501cb1e1211SMatthew G Knepley ierr = PetscSectionSetUp(sectionAdj);CHKERRQ(ierr); 502cb1e1211SMatthew G Knepley if (debug) { 503cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Adjacency Section for Preallocation:\n");CHKERRQ(ierr); 504a8c5bd36SStefano Zampini ierr = PetscSectionView(sectionAdj, NULL);CHKERRQ(ierr); 505cb1e1211SMatthew G Knepley } 506cb1e1211SMatthew G Knepley /* Get adjacent indices */ 507cb1e1211SMatthew G Knepley ierr = PetscSectionGetStorageSize(sectionAdj, &numCols);CHKERRQ(ierr); 508785e854fSJed Brown ierr = PetscMalloc1(numCols, &cols);CHKERRQ(ierr); 509cb1e1211SMatthew G Knepley for (p = pStart; p < pEnd; ++p) { 51076185916SToby Isaac PetscInt numAdj = PETSC_DETERMINE, dof, cdof, off, goff, d, q, anDof, anOff; 511cb1e1211SMatthew G Knepley PetscBool found = PETSC_TRUE; 512cb1e1211SMatthew G Knepley 513cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 514cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 515cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 516cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, p, &goff);CHKERRQ(ierr); 517cb1e1211SMatthew G Knepley for (d = 0; d < dof-cdof; ++d) { 518cb1e1211SMatthew G Knepley PetscInt ldof, rdof; 519cb1e1211SMatthew G Knepley 520cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(leafSectionAdj, off+d, &ldof);CHKERRQ(ierr); 521cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(rootSectionAdj, off+d, &rdof);CHKERRQ(ierr); 522cb1e1211SMatthew G Knepley if (ldof > 0) { 523cb1e1211SMatthew G Knepley /* We do not own this point */ 524cb1e1211SMatthew G Knepley } else if (rdof > 0) { 525cb1e1211SMatthew G Knepley PetscInt aoff, roff; 526cb1e1211SMatthew G Knepley 527cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, goff+d, &aoff);CHKERRQ(ierr); 528cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(rootSectionAdj, off+d, &roff);CHKERRQ(ierr); 529580bdb30SBarry Smith ierr = PetscArraycpy(&cols[aoff], &rootAdj[roff], rdof);CHKERRQ(ierr); 530cb1e1211SMatthew G Knepley } else { 531cb1e1211SMatthew G Knepley found = PETSC_FALSE; 532cb1e1211SMatthew G Knepley } 533cb1e1211SMatthew G Knepley } 534cb1e1211SMatthew G Knepley if (found) continue; 535a9fb6443SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, &numAdj, &tmpAdj);CHKERRQ(ierr); 53676185916SToby Isaac ierr = PetscSectionGetDof(anchorSectionAdj, p, &anDof);CHKERRQ(ierr); 53776185916SToby Isaac ierr = PetscSectionGetOffset(anchorSectionAdj, p, &anOff);CHKERRQ(ierr); 538cb1e1211SMatthew G Knepley for (d = goff; d < goff+dof-cdof; ++d) { 539cb1e1211SMatthew G Knepley PetscInt adof, aoff, i = 0; 540cb1e1211SMatthew G Knepley 541cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, d, &adof);CHKERRQ(ierr); 542cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, d, &aoff);CHKERRQ(ierr); 543cb1e1211SMatthew G Knepley for (q = 0; q < numAdj; ++q) { 544cb1e1211SMatthew G Knepley const PetscInt padj = tmpAdj[q]; 545cb1e1211SMatthew G Knepley PetscInt ndof, ncdof, ngoff, nd; 546cb1e1211SMatthew G Knepley const PetscInt *ncind; 547cb1e1211SMatthew G Knepley 548cb1e1211SMatthew G Knepley /* Adjacent points may not be in the section chart */ 549cb1e1211SMatthew G Knepley if ((padj < pStart) || (padj >= pEnd)) continue; 550cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(section, padj, &ndof);CHKERRQ(ierr); 551cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintDof(section, padj, &ncdof);CHKERRQ(ierr); 552cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstraintIndices(section, padj, &ncind);CHKERRQ(ierr); 553cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionGlobal, padj, &ngoff);CHKERRQ(ierr); 554cb1e1211SMatthew G Knepley for (nd = 0; nd < ndof-ncdof; ++nd, ++i) { 555cb1e1211SMatthew G Knepley cols[aoff+i] = ngoff < 0 ? -(ngoff+1)+nd : ngoff+nd; 556cb1e1211SMatthew G Knepley } 557cb1e1211SMatthew G Knepley } 55876185916SToby Isaac for (q = 0; q < anDof; q++, i++) { 55976185916SToby Isaac cols[aoff+i] = anchorAdj[anOff + q]; 56076185916SToby Isaac } 561cb1e1211SMatthew 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); 562cb1e1211SMatthew G Knepley } 563cb1e1211SMatthew G Knepley } 56476185916SToby Isaac ierr = PetscSectionDestroy(&anchorSectionAdj);CHKERRQ(ierr); 565cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(&leafSectionAdj);CHKERRQ(ierr); 566cb1e1211SMatthew G Knepley ierr = PetscSectionDestroy(&rootSectionAdj);CHKERRQ(ierr); 56776185916SToby Isaac ierr = PetscFree(anchorAdj);CHKERRQ(ierr); 568cb1e1211SMatthew G Knepley ierr = PetscFree(rootAdj);CHKERRQ(ierr); 56970034214SMatthew G. Knepley ierr = PetscFree(tmpAdj);CHKERRQ(ierr); 570cb1e1211SMatthew G Knepley /* Debugging */ 571cb1e1211SMatthew G Knepley if (debug) { 572cb1e1211SMatthew G Knepley IS tmp; 573cb1e1211SMatthew G Knepley ierr = PetscPrintf(comm, "Column indices\n");CHKERRQ(ierr); 574cb1e1211SMatthew G Knepley ierr = ISCreateGeneral(comm, numCols, cols, PETSC_USE_POINTER, &tmp);CHKERRQ(ierr); 575cb1e1211SMatthew G Knepley ierr = ISView(tmp, NULL);CHKERRQ(ierr); 5760a7e1f69SMatthew G. Knepley ierr = ISDestroy(&tmp);CHKERRQ(ierr); 577cb1e1211SMatthew G Knepley } 578a9fb6443SMatthew G. Knepley 579a9fb6443SMatthew G. Knepley *sA = sectionAdj; 580a9fb6443SMatthew G. Knepley *colIdx = cols; 581a9fb6443SMatthew G. Knepley PetscFunctionReturn(0); 582a9fb6443SMatthew G. Knepley } 583a9fb6443SMatthew G. Knepley 58425e402d2SMatthew G. Knepley 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[]) 585a9fb6443SMatthew G. Knepley { 586e101f074SMatthew G. Knepley PetscSection section; 587a9fb6443SMatthew G. Knepley PetscInt rStart, rEnd, r, pStart, pEnd, p; 588a9fb6443SMatthew G. Knepley PetscErrorCode ierr; 589a9fb6443SMatthew G. Knepley 590a9fb6443SMatthew G. Knepley PetscFunctionBegin; 591a9fb6443SMatthew G. Knepley /* This loop needs to change to a loop over points, then field dofs, which means we need to look both sections */ 592cb1e1211SMatthew G Knepley ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 593a9fb6443SMatthew G. Knepley if (rStart%bs || rEnd%bs) SETERRQ3(PetscObjectComm((PetscObject) rLayout), PETSC_ERR_ARG_WRONG, "Invalid layout [%d, %d) for matrix, must be divisible by block size %d", rStart, rEnd, bs); 594a9fb6443SMatthew G. Knepley if (f >= 0 && bs == 1) { 59592fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 596a9fb6443SMatthew G. Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 597a9fb6443SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 598294b7009SMatthew G. Knepley PetscInt rS, rE; 599a9fb6443SMatthew G. Knepley 600088fd00dSMatthew G. Knepley ierr = DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE);CHKERRQ(ierr); 601a9fb6443SMatthew G. Knepley for (r = rS; r < rE; ++r) { 602a9fb6443SMatthew G. Knepley PetscInt numCols, cStart, c; 603a9fb6443SMatthew G. Knepley 604a9fb6443SMatthew G. Knepley ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr); 605a9fb6443SMatthew G. Knepley ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr); 606a9fb6443SMatthew G. Knepley for (c = cStart; c < cStart+numCols; ++c) { 607a9fb6443SMatthew G. Knepley if ((cols[c] >= rStart) && (cols[c] < rEnd)) { 608a9fb6443SMatthew G. Knepley ++dnz[r-rStart]; 609a9fb6443SMatthew G. Knepley if (cols[c] >= r) ++dnzu[r-rStart]; 610a9fb6443SMatthew G. Knepley } else { 611a9fb6443SMatthew G. Knepley ++onz[r-rStart]; 612a9fb6443SMatthew G. Knepley if (cols[c] >= r) ++onzu[r-rStart]; 613a9fb6443SMatthew G. Knepley } 614a9fb6443SMatthew G. Knepley } 615a9fb6443SMatthew G. Knepley } 616a9fb6443SMatthew G. Knepley } 617a9fb6443SMatthew G. Knepley } else { 618cb1e1211SMatthew G Knepley /* Only loop over blocks of rows */ 619cb1e1211SMatthew G Knepley for (r = rStart/bs; r < rEnd/bs; ++r) { 620cb1e1211SMatthew G Knepley const PetscInt row = r*bs; 621cb1e1211SMatthew G Knepley PetscInt numCols, cStart, c; 622cb1e1211SMatthew G Knepley 623cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(sectionAdj, row, &numCols);CHKERRQ(ierr); 624cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(sectionAdj, row, &cStart);CHKERRQ(ierr); 625cb1e1211SMatthew G Knepley for (c = cStart; c < cStart+numCols; ++c) { 626e7bcfa36SMatthew G. Knepley if ((cols[c] >= rStart) && (cols[c] < rEnd)) { 627e7bcfa36SMatthew G. Knepley ++dnz[r-rStart/bs]; 628e7bcfa36SMatthew G. Knepley if (cols[c] >= row) ++dnzu[r-rStart/bs]; 629cb1e1211SMatthew G Knepley } else { 630e7bcfa36SMatthew G. Knepley ++onz[r-rStart/bs]; 631e7bcfa36SMatthew G. Knepley if (cols[c] >= row) ++onzu[r-rStart/bs]; 632cb1e1211SMatthew G Knepley } 633cb1e1211SMatthew G Knepley } 634cb1e1211SMatthew G Knepley } 635a9fb6443SMatthew G. Knepley for (r = 0; r < (rEnd - rStart)/bs; ++r) { 636cb1e1211SMatthew G Knepley dnz[r] /= bs; 637cb1e1211SMatthew G Knepley onz[r] /= bs; 638cb1e1211SMatthew G Knepley dnzu[r] /= bs; 639cb1e1211SMatthew G Knepley onzu[r] /= bs; 640cb1e1211SMatthew G Knepley } 641cb1e1211SMatthew G Knepley } 642a9fb6443SMatthew G. Knepley PetscFunctionReturn(0); 643a9fb6443SMatthew G. Knepley } 644a9fb6443SMatthew G. Knepley 64525e402d2SMatthew G. Knepley static PetscErrorCode DMPlexFillMatrix_Static(DM dm, PetscLayout rLayout, PetscInt bs, PetscInt f, PetscSection sectionAdj, const PetscInt cols[], Mat A) 646a9fb6443SMatthew G. Knepley { 647e101f074SMatthew G. Knepley PetscSection section; 648a9fb6443SMatthew G. Knepley PetscScalar *values; 649a9fb6443SMatthew G. Knepley PetscInt rStart, rEnd, r, pStart, pEnd, p, len, maxRowLen = 0; 650a9fb6443SMatthew G. Knepley PetscErrorCode ierr; 651a9fb6443SMatthew G. Knepley 652a9fb6443SMatthew G. Knepley PetscFunctionBegin; 653a9fb6443SMatthew G. Knepley ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 654a9fb6443SMatthew G. Knepley for (r = rStart; r < rEnd; ++r) { 655a9fb6443SMatthew G. Knepley ierr = PetscSectionGetDof(sectionAdj, r, &len);CHKERRQ(ierr); 656a9fb6443SMatthew G. Knepley maxRowLen = PetscMax(maxRowLen, len); 657a9fb6443SMatthew G. Knepley } 658a9fb6443SMatthew G. Knepley ierr = PetscCalloc1(maxRowLen, &values);CHKERRQ(ierr); 659a9fb6443SMatthew G. Knepley if (f >=0 && bs == 1) { 66092fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 661a9fb6443SMatthew G. Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 662a9fb6443SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 663294b7009SMatthew G. Knepley PetscInt rS, rE; 664a9fb6443SMatthew G. Knepley 665088fd00dSMatthew G. Knepley ierr = DMGetGlobalFieldOffset_Private(dm, p, f, &rS, &rE);CHKERRQ(ierr); 666a9fb6443SMatthew G. Knepley for (r = rS; r < rE; ++r) { 6677e01ee02SMatthew G. Knepley PetscInt numCols, cStart; 668a9fb6443SMatthew G. Knepley 669a9fb6443SMatthew G. Knepley ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr); 670a9fb6443SMatthew G. Knepley ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr); 671a9fb6443SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr); 672a9fb6443SMatthew G. Knepley } 673a9fb6443SMatthew G. Knepley } 674a9fb6443SMatthew G. Knepley } else { 675a9fb6443SMatthew G. Knepley for (r = rStart; r < rEnd; ++r) { 676a9fb6443SMatthew G. Knepley PetscInt numCols, cStart; 677a9fb6443SMatthew G. Knepley 678a9fb6443SMatthew G. Knepley ierr = PetscSectionGetDof(sectionAdj, r, &numCols);CHKERRQ(ierr); 679a9fb6443SMatthew G. Knepley ierr = PetscSectionGetOffset(sectionAdj, r, &cStart);CHKERRQ(ierr); 680a9fb6443SMatthew G. Knepley ierr = MatSetValues(A, 1, &r, numCols, &cols[cStart], values, INSERT_VALUES);CHKERRQ(ierr); 681a9fb6443SMatthew G. Knepley } 682a9fb6443SMatthew G. Knepley } 683a9fb6443SMatthew G. Knepley ierr = PetscFree(values);CHKERRQ(ierr); 684a9fb6443SMatthew G. Knepley PetscFunctionReturn(0); 685a9fb6443SMatthew G. Knepley } 686a9fb6443SMatthew G. Knepley 68725e402d2SMatthew G. Knepley /*@C 68825e402d2SMatthew G. Knepley DMPlexPreallocateOperator - Calculate the matrix nonzero pattern based upon the information in the DM, 68925e402d2SMatthew G. Knepley the PetscDS it contains, and the default PetscSection. 69025e402d2SMatthew G. Knepley 69125e402d2SMatthew G. Knepley Collective 69225e402d2SMatthew G. Knepley 69325e402d2SMatthew G. Knepley Input Arguments: 69425e402d2SMatthew G. Knepley + dm - The DMPlex 69525e402d2SMatthew G. Knepley . bs - The matrix blocksize 69625e402d2SMatthew G. Knepley . dnz - An array to hold the number of nonzeros in the diagonal block 69725e402d2SMatthew G. Knepley . onz - An array to hold the number of nonzeros in the off-diagonal block 69825e402d2SMatthew G. Knepley . dnzu - An array to hold the number of nonzeros in the upper triangle of the diagonal block 69925e402d2SMatthew G. Knepley . onzu - An array to hold the number of nonzeros in the upper triangle of the off-diagonal block 700a2b725a8SWilliam Gropp - fillMatrix - If PETSC_TRUE, fill the matrix with zeros 70125e402d2SMatthew G. Knepley 702*01d2d390SJose E. Roman Output Argument: 70325e402d2SMatthew G. Knepley . A - The preallocated matrix 70425e402d2SMatthew G. Knepley 70525e402d2SMatthew G. Knepley Level: advanced 70625e402d2SMatthew G. Knepley 70725e402d2SMatthew G. Knepley .seealso: DMCreateMatrix() 70825e402d2SMatthew G. Knepley @*/ 709e101f074SMatthew G. Knepley PetscErrorCode DMPlexPreallocateOperator(DM dm, PetscInt bs, PetscInt dnz[], PetscInt onz[], PetscInt dnzu[], PetscInt onzu[], Mat A, PetscBool fillMatrix) 710a9fb6443SMatthew G. Knepley { 711a9fb6443SMatthew G. Knepley MPI_Comm comm; 712a9fb6443SMatthew G. Knepley PetscDS prob; 713a9fb6443SMatthew G. Knepley MatType mtype; 714a9fb6443SMatthew G. Knepley PetscSF sf, sfDof; 715e101f074SMatthew G. Knepley PetscSection section; 716a9fb6443SMatthew G. Knepley PetscInt *remoteOffsets; 717a9fb6443SMatthew G. Knepley PetscSection sectionAdj[4] = {NULL, NULL, NULL, NULL}; 718a9fb6443SMatthew G. Knepley PetscInt *cols[4] = {NULL, NULL, NULL, NULL}; 719a9fb6443SMatthew G. Knepley PetscBool useCone, useClosure; 7207e01ee02SMatthew G. Knepley PetscInt Nf, f, idx, locRows; 721a9fb6443SMatthew G. Knepley PetscLayout rLayout; 722a9fb6443SMatthew G. Knepley PetscBool isSymBlock, isSymSeqBlock, isSymMPIBlock, debug = PETSC_FALSE; 7236cded2eaSMatthew G. Knepley PetscMPIInt size; 724a9fb6443SMatthew G. Knepley PetscErrorCode ierr; 725a9fb6443SMatthew G. Knepley 726a9fb6443SMatthew G. Knepley PetscFunctionBegin; 727a9fb6443SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 728064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(A, MAT_CLASSID, 7); 729064a246eSJacob Faibussowitsch if (dnz) PetscValidPointer(dnz,3); if (onz) PetscValidPointer(onz,4); 730064a246eSJacob Faibussowitsch if (dnzu) PetscValidPointer(dnzu,5); if (onzu) PetscValidPointer(onzu,6); 731a9fb6443SMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 732a9fb6443SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 73392fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 734c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL, "-dm_view_preallocation", &debug, NULL);CHKERRQ(ierr); 735a9fb6443SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 736ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 737a9fb6443SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr); 738a9fb6443SMatthew G. Knepley /* Create dof SF based on point SF */ 739a9fb6443SMatthew G. Knepley if (debug) { 740e101f074SMatthew G. Knepley PetscSection section, sectionGlobal; 741a9fb6443SMatthew G. Knepley PetscSF sf; 742a9fb6443SMatthew G. Knepley 743a9fb6443SMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 74492fd8e1eSJed Brown ierr = DMGetLocalSection(dm, §ion);CHKERRQ(ierr); 745e87a4003SBarry Smith ierr = DMGetGlobalSection(dm, §ionGlobal);CHKERRQ(ierr); 746a9fb6443SMatthew G. Knepley ierr = PetscPrintf(comm, "Input Section for Preallocation:\n");CHKERRQ(ierr); 747a8c5bd36SStefano Zampini ierr = PetscSectionView(section, NULL);CHKERRQ(ierr); 748a9fb6443SMatthew G. Knepley ierr = PetscPrintf(comm, "Input Global Section for Preallocation:\n");CHKERRQ(ierr); 749a8c5bd36SStefano Zampini ierr = PetscSectionView(sectionGlobal, NULL);CHKERRQ(ierr); 7506cded2eaSMatthew G. Knepley if (size > 1) { 751a9fb6443SMatthew G. Knepley ierr = PetscPrintf(comm, "Input SF for Preallocation:\n");CHKERRQ(ierr); 752a9fb6443SMatthew G. Knepley ierr = PetscSFView(sf, NULL);CHKERRQ(ierr); 753a9fb6443SMatthew G. Knepley } 7546cded2eaSMatthew G. Knepley } 755a9fb6443SMatthew G. Knepley ierr = PetscSFCreateRemoteOffsets(sf, section, section, &remoteOffsets);CHKERRQ(ierr); 756a9fb6443SMatthew G. Knepley ierr = PetscSFCreateSectionSF(sf, section, remoteOffsets, section, &sfDof);CHKERRQ(ierr); 7578a713f3bSLawrence Mitchell ierr = PetscFree(remoteOffsets);CHKERRQ(ierr); 7586cded2eaSMatthew G. Knepley if (debug && size > 1) { 759a9fb6443SMatthew G. Knepley ierr = PetscPrintf(comm, "Dof SF for Preallocation:\n");CHKERRQ(ierr); 760a9fb6443SMatthew G. Knepley ierr = PetscSFView(sfDof, NULL);CHKERRQ(ierr); 761a9fb6443SMatthew G. Knepley } 762a9fb6443SMatthew G. Knepley /* Create allocation vectors from adjacency graph */ 763a9fb6443SMatthew G. Knepley ierr = MatGetLocalSize(A, &locRows, NULL);CHKERRQ(ierr); 764a9fb6443SMatthew G. Knepley ierr = PetscLayoutCreate(comm, &rLayout);CHKERRQ(ierr); 765a9fb6443SMatthew G. Knepley ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 766a9fb6443SMatthew G. Knepley ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 767a9fb6443SMatthew G. Knepley ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 768a9fb6443SMatthew G. Knepley /* There are 4 types of adjacency */ 769a9fb6443SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 770a9fb6443SMatthew G. Knepley if (Nf < 1 || bs > 1) { 771b0441da4SMatthew G. Knepley ierr = DMGetBasicAdjacency(dm, &useCone, &useClosure);CHKERRQ(ierr); 772a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 773e101f074SMatthew G. Knepley ierr = DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, §ionAdj[idx], &cols[idx]);CHKERRQ(ierr); 774e101f074SMatthew G. Knepley ierr = DMPlexUpdateAllocation_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu);CHKERRQ(ierr); 775a9fb6443SMatthew G. Knepley } else { 776a9fb6443SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 77734aa8a36SMatthew G. Knepley ierr = DMGetAdjacency(dm, f, &useCone, &useClosure);CHKERRQ(ierr); 778a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 779e101f074SMatthew G. Knepley if (!sectionAdj[idx]) {ierr = DMPlexCreateAdjacencySection_Static(dm, bs, sfDof, useCone, useClosure, PETSC_TRUE, §ionAdj[idx], &cols[idx]);CHKERRQ(ierr);} 780e101f074SMatthew G. Knepley ierr = DMPlexUpdateAllocation_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], dnz, onz, dnzu, onzu);CHKERRQ(ierr); 781a9fb6443SMatthew G. Knepley } 782a9fb6443SMatthew G. Knepley } 783a9fb6443SMatthew G. Knepley ierr = PetscSFDestroy(&sfDof);CHKERRQ(ierr); 784cb1e1211SMatthew G Knepley /* Set matrix pattern */ 785cb1e1211SMatthew G Knepley ierr = MatXAIJSetPreallocation(A, bs, dnz, onz, dnzu, onzu);CHKERRQ(ierr); 786cb1e1211SMatthew G Knepley ierr = MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 78789545effSMatthew G. Knepley /* Check for symmetric storage */ 78889545effSMatthew G. Knepley ierr = MatGetType(A, &mtype);CHKERRQ(ierr); 78989545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATSBAIJ, &isSymBlock);CHKERRQ(ierr); 79089545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATSEQSBAIJ, &isSymSeqBlock);CHKERRQ(ierr); 79189545effSMatthew G. Knepley ierr = PetscStrcmp(mtype, MATMPISBAIJ, &isSymMPIBlock);CHKERRQ(ierr); 79289545effSMatthew G. Knepley if (isSymBlock || isSymSeqBlock || isSymMPIBlock) {ierr = MatSetOption(A, MAT_IGNORE_LOWER_TRIANGULAR, PETSC_TRUE);CHKERRQ(ierr);} 793cb1e1211SMatthew G Knepley /* Fill matrix with zeros */ 794cb1e1211SMatthew G Knepley if (fillMatrix) { 795a9fb6443SMatthew G. Knepley if (Nf < 1 || bs > 1) { 796b0441da4SMatthew G. Knepley ierr = DMGetBasicAdjacency(dm, &useCone, &useClosure);CHKERRQ(ierr); 797a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 798e101f074SMatthew G. Knepley ierr = DMPlexFillMatrix_Static(dm, rLayout, bs, -1, sectionAdj[idx], cols[idx], A);CHKERRQ(ierr); 799a9fb6443SMatthew G. Knepley } else { 800a9fb6443SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 80134aa8a36SMatthew G. Knepley ierr = DMGetAdjacency(dm, f, &useCone, &useClosure);CHKERRQ(ierr); 802a9fb6443SMatthew G. Knepley idx = (useCone ? 1 : 0) + (useClosure ? 2 : 0); 803e101f074SMatthew G. Knepley ierr = DMPlexFillMatrix_Static(dm, rLayout, bs, f, sectionAdj[idx], cols[idx], A);CHKERRQ(ierr); 804cb1e1211SMatthew G Knepley } 805cb1e1211SMatthew G Knepley } 806cb1e1211SMatthew G Knepley ierr = MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 807cb1e1211SMatthew G Knepley ierr = MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 808cb1e1211SMatthew G Knepley } 809a9fb6443SMatthew G. Knepley ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 810a9fb6443SMatthew G. Knepley for (idx = 0; idx < 4; ++idx) {ierr = PetscSectionDestroy(§ionAdj[idx]);CHKERRQ(ierr); ierr = PetscFree(cols[idx]);CHKERRQ(ierr);} 811a72f3261SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Preallocate,dm,0,0,0);CHKERRQ(ierr); 812cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 813cb1e1211SMatthew G Knepley } 814cb1e1211SMatthew G Knepley 815cb1e1211SMatthew G Knepley #if 0 816cb1e1211SMatthew 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) 817cb1e1211SMatthew G Knepley { 818cb1e1211SMatthew G Knepley PetscInt *tmpClosure,*tmpAdj,*visits; 819cb1e1211SMatthew G Knepley PetscInt c,cStart,cEnd,pStart,pEnd; 820cb1e1211SMatthew G Knepley PetscErrorCode ierr; 821cb1e1211SMatthew G Knepley 822cb1e1211SMatthew G Knepley PetscFunctionBegin; 823c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 824cb1e1211SMatthew G Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 825cb1e1211SMatthew G Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 826cb1e1211SMatthew G Knepley 827e7b80042SMatthew G. Knepley maxClosureSize = 2*PetscMax(PetscPowInt(mesh->maxConeSize,depth+1),PetscPowInt(mesh->maxSupportSize,depth+1)); 828cb1e1211SMatthew G Knepley 829cb1e1211SMatthew G Knepley ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 830cb1e1211SMatthew G Knepley npoints = pEnd - pStart; 831cb1e1211SMatthew G Knepley 832dcca6d9dSJed Brown ierr = PetscMalloc3(maxClosureSize,&tmpClosure,npoints,&lvisits,npoints,&visits);CHKERRQ(ierr); 833580bdb30SBarry Smith ierr = PetscArrayzero(lvisits,pEnd-pStart);CHKERRQ(ierr); 834580bdb30SBarry Smith ierr = PetscArrayzero(visits,pEnd-pStart);CHKERRQ(ierr); 835cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 836cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 837cb1e1211SMatthew G Knepley PetscInt *support = tmpClosure; 838cb1e1211SMatthew G Knepley ierr = DMPlexGetTransitiveClosure(dm, c, PETSC_FALSE, &supportSize, (PetscInt**)&support);CHKERRQ(ierr); 839cb1e1211SMatthew G Knepley for (p=0; p<supportSize; p++) lvisits[support[p]]++; 840cb1e1211SMatthew G Knepley } 841cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr); 842cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd (sf,MPIU_INT,lvisits,visits,MPI_SUM);CHKERRQ(ierr); 843ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(sf,MPIU_INT,visits,lvisits,MPI_REPLACE);CHKERRQ(ierr); 844cb1e1211SMatthew G Knepley ierr = PetscSFBcastEnd (sf,MPIU_INT,visits,lvisits);CHKERRQ(ierr); 845cb1e1211SMatthew G Knepley 846dec1416fSJunchao Zhang ierr = PetscSFGetRootRanks();CHKERRQ(ierr); 847cb1e1211SMatthew G Knepley 848dcca6d9dSJed Brown ierr = PetscMalloc2(maxClosureSize*maxClosureSize,&cellmat,npoints,&owner);CHKERRQ(ierr); 849cb1e1211SMatthew G Knepley for (c=cStart; c<cEnd; c++) { 850580bdb30SBarry Smith ierr = PetscArrayzero(cellmat,maxClosureSize*maxClosureSize);CHKERRQ(ierr); 851cb1e1211SMatthew G Knepley /* 852cb1e1211SMatthew G Knepley Depth-first walk of transitive closure. 853cb1e1211SMatthew 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. 854cb1e1211SMatthew G Knepley This contribution is added to dnz if owning ranks of p and q match, to onz otherwise. 855cb1e1211SMatthew G Knepley */ 856cb1e1211SMatthew G Knepley } 857cb1e1211SMatthew G Knepley 858cb1e1211SMatthew G Knepley ierr = PetscSFReduceBegin(sf,MPIU_INT,ldnz,dnz,MPI_SUM);CHKERRQ(ierr); 859cb1e1211SMatthew G Knepley ierr = PetscSFReduceEnd (sf,MPIU_INT,lonz,onz,MPI_SUM);CHKERRQ(ierr); 860cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 861cb1e1211SMatthew G Knepley } 862cb1e1211SMatthew G Knepley #endif 863