Lines Matching +full:petsc +full:- +full:dist

1 #include <petsc/private/dmpleximpl.h>  /*I      "petscdmplex.h"   I*/
2 #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/
3 #include <petsc/private/partitionerimpl.h>
6 DMPlexSetAdjacencyUser - Define adjacency in the mesh using a user-provided callback
9 + dm - The DM object
10 . user - The user callback, may be `NULL` (to clear the callback)
11 - ctx - context for callback evaluation, may be `NULL`
24 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexSetAdjacencyUser()
28 mesh->useradjacency = user; in DMPlexSetAdjacencyUser()
29 mesh->useradjacencyctx = ctx; in DMPlexSetAdjacencyUser()
34 DMPlexGetAdjacencyUser - get the user-defined adjacency callback
37 . dm - The `DM` object
40 + user - The callback
41 - ctx - context for callback evaluation
49 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexGetAdjacencyUser()
53 if (user) *user = mesh->useradjacency; in DMPlexGetAdjacencyUser()
54 if (ctx) *ctx = mesh->useradjacencyctx; in DMPlexGetAdjacencyUser()
59 DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
62 + dm - The `DM` object
63 - useAnchors - Flag to use the constraints. If PETSC_TRUE, then constrained points are omitted fro…
71 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexSetAdjacencyUseAnchors()
75 mesh->useAnchors = useAnchors; in DMPlexSetAdjacencyUseAnchors()
80 …DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constra…
83 . dm - The `DM` object
86 . useAnchors - Flag to use the closure. If PETSC_TRUE, then constrained points are omitted from DM…
94 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexGetAdjacencyUseAnchors()
99 *useAnchors = mesh->useAnchors; in DMPlexGetAdjacencyUseAnchors()
112 const PetscInt point = !c ? p : cone[c - 1]; in DMPlexGetAdjacency_Cone_Internal()
138 const PetscInt point = !s ? p : support[s - 1]; in DMPlexGetAdjacency_Support_Internal()
200 depth = PetscMax(depth, -depth); in DMPlexGetMaxAdjacencySize_Internal()
206 = \sum^d_{i=0} (maxS*maxC)^i - 1 in DMPlexGetMaxAdjacencySize_Internal()
207 = (maxS*maxC)^{d+1} - 1 / (maxS*maxC - 1) - 1 in DMPlexGetMaxAdjacencySize_Internal()
209 …supp[d](cell) + supp[d-1](maxC[d] faces) + supp[1](maxC[d]*maxC[d-1] edges) + supp[0](maxC[d]*maxC… in DMPlexGetMaxAdjacencySize_Internal()
210 = 0 + maxS[d-1]*maxC[d] + maxS[1]*maxC[d]*maxC[d-1] + maxS[0]*maxC[d]*maxC[d-1]*maxC[d-2] in DMPlexGetMaxAdjacencySize_Internal()
213 if ((depth == 3 && maxP > 200) || (depth == 2 && maxP > 580)) asiz = pEnd - pStart; in DMPlexGetMaxAdjacencySize_Internal()
214 else asiz = (maxP > 1) ? ((PetscPowInt(maxP, depth + 1) - 1) / (maxP - 1)) : depth + 1; in DMPlexGetMaxAdjacencySize_Internal()
216 *max_adjacency_size = PetscMin(asiz, pEnd - pStart); in DMPlexGetMaxAdjacencySize_Internal()
222 // + adjSize - Number of adjacent points
223 // - adj - Array of the adjacent points
227 PetscInt aStart = -1, aEnd = -1; in DMPlexGetAdjacency_Internal()
232 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexGetAdjacency_Internal()
248 if (mesh->useradjacency) { in DMPlexGetAdjacency_Internal()
249 PetscCall((*mesh->useradjacency)(dm, p, adjSize, *adj, mesh->useradjacencyctx)); in DMPlexGetAdjacency_Internal()
272 for (j = i + 1; j < numAdj; j++) orig[j - 1] = orig[j]; in DMPlexGetAdjacency_Internal()
273 origSize--; in DMPlexGetAdjacency_Internal()
274 numAdj--; in DMPlexGetAdjacency_Internal()
293 DMPlexGetAdjacency - Return all points adjacent to the given point
296 + dm - The `DM` object
297 - p - The point
300 + adjSize - The maximum size of `adj` if it is non-`NULL`, or `PETSC_DETERMINE`;
302 - adj - Either `NULL` so that the array is allocated, or an existing array with size `adjSize`;
327 DMPlexCreateTwoSidedProcessSF - Create an `PetscSF` which just has process connectivity
332 + dm - The `DM`
333 . sfPoint - The `PetscSF` which encodes point connectivity
334 . rootRankSection - to be documented
335 . rootRanks - to be documented
336 . leafRankSection - to be documented
337 - leafRanks - to be documented
340 + processRanks - A list of process neighbors, or `NULL`
341 - sfProcess - An `PetscSF` encoding the two-sided process connectivity, or `NULL`
368 /* Compute root-to-leaf process connectivity */ in DMPlexCreateTwoSidedProcessSF()
379 /* Compute leaf-to-neighbor process connectivity */ in DMPlexCreateTwoSidedProcessSF()
390 /* Compute leaf-to-root process connectivity */ in DMPlexCreateTwoSidedProcessSF()
414 PetscCall(PetscObjectSetName((PetscObject)*sfProcess, "Two-Sided Process SF")); in DMPlexCreateTwoSidedProcessSF()
422 …DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-s…
427 . dm - The `DM`
430 + rootSection - The number of leaves for a given root point
431 . rootrank - The rank of each edge into the root point
432 . leafSection - The number of processes sharing a given leaf point
433 - leafrank - The rank of each process sharing a leaf point
458 …for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(rootSection, p, rootdegree[p - pStart… in DMPlexDistributeOwnership()
462 PetscCall(PetscMalloc1(pEnd - pStart, &myrank)); in DMPlexDistributeOwnership()
464 for (p = 0; p < pEnd - pStart; ++p) myrank[p] = rank; in DMPlexDistributeOwnership()
476 …DMPlexCreateOverlapLabel - Compute a label indicating what overlap points should be sent to new pr…
481 + dm - The `DM`
482 . levels - Number of overlap levels
483 . rootSection - The number of leaves for a given root point
484 . rootrank - The rank of each edge into the root point
485 . leafSection - The number of processes sharing a given leaf point
486 - leafrank - The rank of each process sharing a leaf point
489 . ovLabel - `DMLabel` containing remote overlap contributions as point/rank pairings
576 …PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-overlap_vie… in DMPlexCreateOverlapLabel()
634 …DMPlexCreateOverlapLabelFromLabels - Compute a label indicating what overlap points should be sent…
639 + dm - The `DM`
640 . numLabels - The number of labels to draw candidate points from
641 . label - An array of labels containing candidate points
642 . value - An array of label values marking the candidate points
643 . numExLabels - The number of labels to use for exclusion
644 . exLabel - An array of labels indicating points to be excluded, or `NULL`
645 . exValue - An array of label values to be excluded, or `NULL`
646 . rootSection - The number of leaves for a given root point
647 . rootrank - The rank of each edge into the root point
648 . leafSection - The number of processes sharing a given leaf point
649 - leafrank - The rank of each process sharing a leaf point
652 . ovLabel - `DMLabel` containing remote overlap contributions as point/rank pairings
704 else loc = (p >= 0 && p < nleaves) ? p : -1; in DMPlexCreateOverlapLabelFromLabels()
737 …PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-overlap_vie… in DMPlexCreateOverlapLabelFromLabels()
758 …DMPlexCreateOverlapMigrationSF - Create a `PetscSF` describing the new mesh distribution to make t…
763 + dm - The `DM`
764 - overlapSF - The `PetscSF` mapping ghost points in overlap to owner points on other processes
767 . migrationSF - A `PetscSF` that maps original points in old locations to points in new locations
799 for (p = 0; p < nleaves; p++) remoteDepths[p] = -1; in DMPlexCreateOverlapMigrationSF()
810 for (d = dim - 2; d > 0; d--) depthShift[d] += depthRecv[d + 1]; in DMPlexCreateOverlapMigrationSF()
818 newLeaves = pEnd - pStart + nleaves; in DMPlexCreateOverlapMigrationSF()
861 …PetscCall(PetscSFSetGraph(*migrationSF, pEnd - pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremo… in DMPlexCreateOverlapMigrationSF()
868 DMPlexStratifyMigrationSF - Rearrange the leaves of a migration sf for stratification.
871 + dm - The DM
872 - sf - A star forest with non-ordered leaves, usually defining a DM point migration
875 . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratifi…
920 remoteDepths[p].index = -1; in DMPlexStratifyMigrationSF()
921 remoteDepths[p].rank = -1; in DMPlexStratifyMigrationSF()
937 /* Cells (depth), Vertices (0), Faces (depth-1), Edges (1) in DMPlexStratifyMigrationSF()
942 depths[2] = depth - 1; in DMPlexStratifyMigrationSF()
946 dims[2] = dim - 1; in DMPlexStratifyMigrationSF()
987 …DMPlexDistributeField - Distribute field data to match a given `PetscSF`, usually the `PetscSF` fr…
992 + dm - The `DMPLEX` object
993 . pointSF - The `PetscSF` describing the communication pattern
994 . originalSection - The `PetscSection` for existing data layout
995 - originalVec - The existing data in a local vector
998 + newSection - The `PetscSF` describing the new data layout
999 - newVec - The new data in a local vector
1017 PetscCall(VecSetType(newVec, dm->vectype)); in DMPlexDistributeField()
1033 …DMPlexDistributeFieldIS - Distribute field data to match a given `PetscSF`, usually the `PetscSF` …
1038 + dm - The `DMPLEX` object
1039 . pointSF - The `PetscSF` describing the communication pattern
1040 . originalSection - The `PetscSection` for existing data layout
1041 - originalIS - The existing data
1044 + newSection - The `PetscSF` describing the new data layout
1045 - newIS - The new data
1069 …DMPlexDistributeData - Distribute field data to match a given `PetscSF`, usually the `PetscSF` fro…
1074 + dm - The `DMPLEX` object
1075 . pointSF - The `PetscSF` describing the communication pattern
1076 . originalSection - The `PetscSection` for existing data layout
1077 . datatype - The type of data
1078 - originalData - The existing data
1081 + newSection - The `PetscSection` describing the new data layout
1082 - newData - The new data
1087 This is simply a wrapper around `PetscSectionMigrateData()`, but includes DM-specific logging.
1102 DM_Plex *pmesh = (DM_Plex *)dmParallel->data; in DMPlexDistributeCones()
1149 …PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-cones_view"… in DMPlexDistributeCones()
1167 PetscCall(PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd)); in DMPlexDistributeCones()
1168 PetscCall(PetscSectionSetChart(pmesh->supportSection, pStart, pEnd)); in DMPlexDistributeCones()
1247 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexDistributeLabels()
1253 PetscObjectState depthState = -1; in DMPlexDistributeLabels()
1267 lsendDepth = mesh->depthState != depthState ? PETSC_TRUE : PETSC_FALSE; in DMPlexDistributeLabels()
1270 PetscCall(DMPlexGetDepthLabel(dmParallel, &dmParallel->depthLabel)); in DMPlexDistributeLabels()
1271 PetscCall(DMRemoveLabelBySelf(dmParallel, &dmParallel->depthLabel, PETSC_FALSE)); in DMPlexDistributeLabels()
1328 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexDistributeSetupTree()
1329 DM_Plex *pmesh = (DM_Plex *)dmParallel->data; in DMPlexDistributeSetupTree()
1385 …PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-parents_vie… in DMPlexDistributeSetupTree()
1398 pmesh->useAnchors = mesh->useAnchors; in DMPlexDistributeSetupTree()
1403 …DMPlexSetPartitionBalance - Should distribution of the `DM` attempt to balance the shared point pa…
1406 + dm - The `DMPLEX` object
1407 - flg - Balance the partition?
1415 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexSetPartitionBalance()
1418 mesh->partitionBalance = flg; in DMPlexSetPartitionBalance()
1423 …DMPlexGetPartitionBalance - Does distribution of the `DM` attempt to balance the shared point part…
1426 . dm - The `DMPLEX` object
1429 . flg - Balance the partition?
1437 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexGetPartitionBalance()
1442 *flg = mesh->partitionBalance; in DMPlexGetPartitionBalance()
1471 DMPlexCreatePointSF - Build a point `PetscSF` from an `PetscSF` describing a point migration
1474 + dm - The source `DMPLEX` object
1475 . migrationSF - The star forest that describes the parallel point remapping
1476 - ownership - Flag causing a vote to determine point ownership
1479 . pointSF - The star forest describing the point overlap in the remapped `DM`
1548 rootVote[p].vote = -3; in DMPlexCreatePointSF()
1549 rootVote[p].rank = -3; in DMPlexCreatePointSF()
1550 rootVote[p].index = -3; in DMPlexCreatePointSF()
1567 rootNodes[p].index = -1; in DMPlexCreatePointSF()
1585 /* Note that pointLocal is automatically sorted as it is sublist of 0, ..., nleaves-1 */ in DMPlexCreatePointSF()
1604 DMPlexMigrate - Migrates internal `DM` data over the supplied star forest
1609 + dm - The source `DMPLEX` object
1610 - sf - The star forest communication context describing the migration pattern
1613 . targetDM - The target `DMPLEX` object
1636 /* Check for a one-to-all distribution pattern */ in DMPlexMigrate()
1651 if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n] + 1); in DMPlexMigrate()
1653 /* Derive the new local-to-global mapping from the old one */ in DMPlexMigrate()
1662 /* One-to-all distribution pattern: We can derive LToG from SF */ in DMPlexMigrate()
1666 …PetscCall(ISLocalToGlobalMappingViewFromOptions(ltogMigration, (PetscObject)dm, "-partition_view")… in DMPlexMigrate()
1679 DMPlexRemapMigrationSF - Rewrite the distribution SF to account for overlap
1684 + sfOverlap - The `PetscSF` object just for the overlap
1685 - sfMigration - The original distribution `PetscSF` object
1688 . sfMigrationNew - A rewritten `PetscSF` object that incorporates the overlap
1721 /* For DG-like methods, the code below is equivalent (but faster) than calling
1735 PetscCall(PetscMalloc1((2 * (pEnd - pStart)), &closure)); in DMPlexCreateClosureIndex_CELL()
1742 …PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 2 * (pEnd - pStart), closure, PETSC_OWN_POINTER, &clPoi… in DMPlexCreateClosureIndex_CELL()
1761 printHeader = part->printHeader; in DMPlexDistribute_Multistage()
1767 PetscInt ovl = (l < nl - 1) ? 0 : overlap; in DMPlexDistribute_Multistage()
1772 PetscCall(DMViewFromOptions(dm, (PetscObject)part, "-petscpartitioner_multistage_dm_view")); in DMPlexDistribute_Multistage()
1778 part->printHeader = PETSC_FALSE; in DMPlexDistribute_Multistage()
1781 if (part->usevwgt && dm->localSection && l != nl - 1) { in DMPlexDistribute_Multistage()
1804 part->printHeader = printHeader; in DMPlexDistribute_Multistage()
1809 DMPlexDistribute - Distributes the mesh and any associated sections.
1814 + dm - The original `DMPLEX` object
1815 - overlap - The overlap of partitions, 0 is the default
1818 + sf - The `PetscSF` used for point distribution, or `NULL` if not needed
1819 - dmParallel - The distributed `DMPLEX` object
1919 …PetscCall(PetscOptionsHasName(((PetscObject)dm)->options, ((PetscObject)dm)->prefix, "-partition_v… in DMPlexDistribute()
1925 /* Create non-overlapping parallel DM and migrate internal data */ in DMPlexDistribute()
1953 /* Re-map the migration SF to establish the full migration pattern */ in DMPlexDistribute()
1967 if (dm->localSection) { in DMPlexDistribute()
1970 PetscCall(PetscSFDistributeSection(sfMigration, dm->localSection, NULL, psection)); in DMPlexDistribute()
1974 if (dm->useNatural) { in DMPlexDistribute()
1983 if (dm->sfNatural) { in DMPlexDistribute()
1984 …cCall(DMPlexMigrateGlobalToNaturalSF(dm, *dmParallel, dm->sfNatural, sfMigration, &(*dmParallel)->… in DMPlexDistribute()
1986 …PetscCall(DMPlexCreateGlobalToNaturalSF(*dmParallel, section, sfMigration, &(*dmParallel)->sfNatur… in DMPlexDistribute()
1989 if (dm->sfMigration) { in DMPlexDistribute()
1992 PetscCall(PetscSFCompose(dm->sfMigration, sfMigration, &naturalPointSF)); in DMPlexDistribute()
2008 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexDistributeOverlap_Internal()
2059->numOvLabels) PetscCall(DMPlexCreateOverlapLabelFromLabels(dm, mesh->numOvLabels, mesh->ovLabels,… in DMPlexDistributeOverlap_Internal()
2111 DMPlexDistributeOverlap - Add partition overlap to a distributed non-overlapping `DM`.
2116 + dm - The non-overlapping distributed `DMPLEX` object
2117 - overlap - The overlap of partitions (the same on all ranks)
2120 + sf - The `PetscSF` used for point distribution, or pass `NULL` if not needed
2121 - dmOverlap - The overlapping distributed `DMPLEX` object
2124 + -dm_plex_overlap_labels <name1,name2,...> - List of overlap label names
2125 . -dm_plex_overlap_values <int1,int2,...> - List of overlap label values
2126 . -dm_plex_overlap_exclude_label <name> - Label used to exclude points from overlap
2127 - -dm_plex_overlap_exclude_value <int> - Label value used to exclude points from overlap
2152 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexGetOverlap_Plex()
2155 *overlap = mesh->overlap; in DMPlexGetOverlap_Plex()
2166 mesh = (DM_Plex *)dm->data; in DMPlexSetOverlap_Plex()
2169 meshSrc = (DM_Plex *)dmSrc->data; in DMPlexSetOverlap_Plex()
2170 mesh->overlap = meshSrc->overlap; in DMPlexSetOverlap_Plex()
2172 mesh->overlap = 0; in DMPlexSetOverlap_Plex()
2174 mesh->overlap += overlap; in DMPlexSetOverlap_Plex()
2179 DMPlexGetOverlap - Get the width of the cell overlap
2184 . dm - The `DM`
2187 . overlap - the width of the cell overlap
2203 DMPlexSetOverlap - Set the width of the cell overlap
2208 + dm - The `DM`
2209 . dmSrc - The `DM` that produced this one, or `NULL`
2210 - overlap - the width of the cell overlap
2228 PetscErrorCode DMPlexDistributeSetDefault_Plex(DM dm, PetscBool dist) in DMPlexDistributeSetDefault_Plex() argument
2230 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexDistributeSetDefault_Plex()
2233 mesh->distDefault = dist; in DMPlexDistributeSetDefault_Plex()
2238 DMPlexDistributeSetDefault - Set flag indicating whether the `DM` should be distributed by default
2243 + dm - The `DM`
2244 - dist - Flag for distribution
2250 PetscErrorCode DMPlexDistributeSetDefault(DM dm, PetscBool dist) in DMPlexDistributeSetDefault() argument
2254 PetscValidLogicalCollectiveBool(dm, dist, 2); in DMPlexDistributeSetDefault()
2255 PetscTryMethod(dm, "DMPlexDistributeSetDefault_C", (DM, PetscBool), (dm, dist)); in DMPlexDistributeSetDefault()
2259 PetscErrorCode DMPlexDistributeGetDefault_Plex(DM dm, PetscBool *dist) in DMPlexDistributeGetDefault_Plex() argument
2261 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexDistributeGetDefault_Plex()
2264 *dist = mesh->distDefault; in DMPlexDistributeGetDefault_Plex()
2269 DMPlexDistributeGetDefault - Get flag indicating whether the `DM` should be distributed by default
2274 . dm - The `DM`
2277 . dist - Flag for distribution
2283 PetscErrorCode DMPlexDistributeGetDefault(DM dm, PetscBool *dist) in DMPlexDistributeGetDefault() argument
2287 PetscAssertPointer(dist, 2); in DMPlexDistributeGetDefault()
2288 PetscUseMethod(dm, "DMPlexDistributeGetDefault_C", (DM, PetscBool *), (dm, dist)); in DMPlexDistributeGetDefault()
2293 DMPlexGetGatherDM - Get a copy of the `DMPLEX` that gathers all points on the
2299 . dm - the original `DMPLEX` object
2302 + sf - the `PetscSF` used for point distribution (optional)
2303 - gatherMesh - the gathered `DM` object, or `NULL`
2337 DMPlexGetRedundantDM - Get a copy of the `DMPLEX` that is completely copied on each process.
2342 . dm - the original `DMPLEX` object
2345 + sf - the `PetscSF` used for point distribution (optional)
2346 - redundantMesh - the redundant `DM` object, or `NULL`
2357 PetscInt numPoints = -1; in DMPlexGetRedundantDM()
2378 numPoints = pEnd - pStart; in DMPlexGetRedundantDM()
2386 …PetscCall(PetscSFSetGraph(migrationSF, pEnd - pStart, numPoints, NULL, PETSC_OWN_POINTER, points, … in DMPlexGetRedundantDM()
2413 …DMPlexIsDistributed - Find out whether this `DM` is distributed, i.e. more than one rank owns some…
2418 . dm - The `DM` object
2421 . distributed - Flag whether the `DM` is distributed
2448 count = (pEnd - pStart) > 0 ? 1 : 0; in DMPlexIsDistributed()
2455 DMPlexDistributionSetName - Set the name of the specific parallel distribution
2458 + dm - The `DM`
2459 - name - The name of the specific parallel distribution
2473 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexDistributionSetName()
2478 PetscCall(PetscFree(mesh->distributionName)); in DMPlexDistributionSetName()
2479 PetscCall(PetscStrallocpy(name, &mesh->distributionName)); in DMPlexDistributionSetName()
2484 DMPlexDistributionGetName - Retrieve the name of the specific parallel distribution
2487 . dm - The `DM`
2490 . name - The name of the specific parallel distribution
2504 DM_Plex *mesh = (DM_Plex *)dm->data; in DMPlexDistributionGetName()
2509 *name = mesh->distributionName; in DMPlexDistributionGetName()