1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/ 2abe9303eSLisandro Dalcin 3abe9303eSLisandro Dalcin typedef struct { 4abe9303eSLisandro Dalcin PetscSection section; /* Sizes for each partition */ 5abe9303eSLisandro Dalcin IS partition; /* Points in each partition */ 6abe9303eSLisandro Dalcin PetscBool random; /* Flag for a random partition */ 7abe9303eSLisandro Dalcin } PetscPartitioner_Shell; 8abe9303eSLisandro Dalcin 9d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerReset_Shell(PetscPartitioner part) 10d71ae5a4SJacob Faibussowitsch { 11abe9303eSLisandro Dalcin PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data; 12abe9303eSLisandro Dalcin 13abe9303eSLisandro Dalcin PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&p->section)); 159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&p->partition)); 163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17abe9303eSLisandro Dalcin } 18abe9303eSLisandro Dalcin 19d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part) 20d71ae5a4SJacob Faibussowitsch { 21abe9303eSLisandro Dalcin PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(PetscPartitionerReset_Shell(part)); 239566063dSJacob Faibussowitsch PetscCall(PetscFree(part->data)); 243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25abe9303eSLisandro Dalcin } 26abe9303eSLisandro Dalcin 27d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Shell_ASCII(PetscPartitioner part, PetscViewer viewer) 28d71ae5a4SJacob Faibussowitsch { 29abe9303eSLisandro Dalcin PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data; 30abe9303eSLisandro Dalcin 31abe9303eSLisandro Dalcin PetscFunctionBegin; 32abe9303eSLisandro Dalcin if (p->random) { 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "using random partition\n")); 359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 36abe9303eSLisandro Dalcin } 373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38abe9303eSLisandro Dalcin } 39abe9303eSLisandro Dalcin 40d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer) 41d71ae5a4SJacob Faibussowitsch { 42abe9303eSLisandro Dalcin PetscBool iascii; 43abe9303eSLisandro Dalcin 44abe9303eSLisandro Dalcin PetscFunctionBegin; 45abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 46abe9303eSLisandro Dalcin PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 489566063dSJacob Faibussowitsch if (iascii) PetscCall(PetscPartitionerView_Shell_ASCII(part, viewer)); 493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50abe9303eSLisandro Dalcin } 51abe9303eSLisandro Dalcin 52d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscPartitioner part, PetscOptionItems *PetscOptionsObject) 53d71ae5a4SJacob Faibussowitsch { 54d7cc930eSLisandro Dalcin PetscBool random = PETSC_FALSE, set; 55abe9303eSLisandro Dalcin 56abe9303eSLisandro Dalcin PetscFunctionBegin; 57d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "PetscPartitioner Shell Options"); 589566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellGetRandom(part, &random)); 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &random, &set)); 609566063dSJacob Faibussowitsch if (set) PetscCall(PetscPartitionerShellSetRandom(part, random)); 61d0609cedSBarry Smith PetscOptionsHeadEnd(); 623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 63abe9303eSLisandro Dalcin } 64abe9303eSLisandro Dalcin 65d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) 66d71ae5a4SJacob Faibussowitsch { 67abe9303eSLisandro Dalcin PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data; 68abe9303eSLisandro Dalcin PetscInt np; 69abe9303eSLisandro Dalcin 70abe9303eSLisandro Dalcin PetscFunctionBegin; 71abe9303eSLisandro Dalcin if (p->random) { 72abe9303eSLisandro Dalcin PetscRandom r; 73abe9303eSLisandro Dalcin PetscInt *sizes, *points, v, p; 74abe9303eSLisandro Dalcin PetscMPIInt rank; 75abe9303eSLisandro Dalcin 769566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)part), &rank)); 779566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r)); 789566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(r, 0.0, (PetscScalar)nparts)); 799566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 809566063dSJacob Faibussowitsch PetscCall(PetscCalloc2(nparts, &sizes, numVertices, &points)); 81ad540459SPierre Jolivet for (v = 0; v < numVertices; ++v) points[v] = v; 82ad540459SPierre Jolivet for (p = 0; p < nparts; ++p) sizes[p] = numVertices / nparts + (PetscInt)(p < numVertices % nparts); 83abe9303eSLisandro Dalcin for (v = numVertices - 1; v > 0; --v) { 84abe9303eSLisandro Dalcin PetscReal val; 85abe9303eSLisandro Dalcin PetscInt w, tmp; 86abe9303eSLisandro Dalcin 879566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(r, 0.0, (PetscScalar)(v + 1))); 889566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValueReal(r, &val)); 89abe9303eSLisandro Dalcin w = PetscFloorReal(val); 90abe9303eSLisandro Dalcin tmp = points[v]; 91abe9303eSLisandro Dalcin points[v] = points[w]; 92abe9303eSLisandro Dalcin points[w] = tmp; 93abe9303eSLisandro Dalcin } 949566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 959566063dSJacob Faibussowitsch PetscCall(PetscPartitionerShellSetPartition(part, nparts, sizes, points)); 969566063dSJacob Faibussowitsch PetscCall(PetscFree2(sizes, points)); 97abe9303eSLisandro Dalcin } 9828b400f6SJacob Faibussowitsch PetscCheck(p->section, PetscObjectComm((PetscObject)part), PETSC_ERR_ARG_WRONG, "Shell partitioner information not provided. Please call PetscPartitionerShellSetPartition()"); 999566063dSJacob Faibussowitsch PetscCall(PetscSectionGetChart(p->section, NULL, &np)); 10063a3b9bcSJacob Faibussowitsch PetscCheck(nparts == np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %" PetscInt_FMT " != configured partitions %" PetscInt_FMT, nparts, np); 1019566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(p->partition, &np)); 10263a3b9bcSJacob Faibussowitsch PetscCheck(numVertices == np, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %" PetscInt_FMT " != configured vertices %" PetscInt_FMT, numVertices, np); 1039566063dSJacob Faibussowitsch PetscCall(PetscSectionCopy(p->section, partSection)); 104abe9303eSLisandro Dalcin *partition = p->partition; 1059566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)p->partition)); 1063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107abe9303eSLisandro Dalcin } 108abe9303eSLisandro Dalcin 109d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part) 110d71ae5a4SJacob Faibussowitsch { 111abe9303eSLisandro Dalcin PetscFunctionBegin; 112abe9303eSLisandro Dalcin part->noGraph = PETSC_TRUE; /* PetscPartitionerShell cannot overload the partition call, so it is safe for now */ 113abe9303eSLisandro Dalcin part->ops->view = PetscPartitionerView_Shell; 114abe9303eSLisandro Dalcin part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell; 115abe9303eSLisandro Dalcin part->ops->reset = PetscPartitionerReset_Shell; 116abe9303eSLisandro Dalcin part->ops->destroy = PetscPartitionerDestroy_Shell; 117abe9303eSLisandro Dalcin part->ops->partition = PetscPartitionerPartition_Shell; 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119abe9303eSLisandro Dalcin } 120abe9303eSLisandro Dalcin 121abe9303eSLisandro Dalcin /*MC 122abe9303eSLisandro Dalcin PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object 123abe9303eSLisandro Dalcin 124abe9303eSLisandro Dalcin Level: intermediate 125abe9303eSLisandro Dalcin 126abe9303eSLisandro Dalcin Options Database Keys: 127abe9303eSLisandro Dalcin . -petscpartitioner_shell_random - Use a random partition 128abe9303eSLisandro Dalcin 129db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()` 130abe9303eSLisandro Dalcin M*/ 131abe9303eSLisandro Dalcin 132d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part) 133d71ae5a4SJacob Faibussowitsch { 134abe9303eSLisandro Dalcin PetscPartitioner_Shell *p; 135abe9303eSLisandro Dalcin 136abe9303eSLisandro Dalcin PetscFunctionBegin; 137abe9303eSLisandro Dalcin PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1); 1384dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&p)); 139abe9303eSLisandro Dalcin part->data = p; 140abe9303eSLisandro Dalcin 1419566063dSJacob Faibussowitsch PetscCall(PetscPartitionerInitialize_Shell(part)); 142abe9303eSLisandro Dalcin p->random = PETSC_FALSE; 1433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 144abe9303eSLisandro Dalcin } 145abe9303eSLisandro Dalcin 146abe9303eSLisandro Dalcin /*@C 14735cb6cd3SPierre Jolivet PetscPartitionerShellSetPartition - Set an artificial partition for a mesh 148abe9303eSLisandro Dalcin 14920f4b53cSBarry Smith Collective 150abe9303eSLisandro Dalcin 151abe9303eSLisandro Dalcin Input Parameters: 15220f4b53cSBarry Smith + part - The `PetscPartitioner` 153abe9303eSLisandro Dalcin . size - The number of partitions 15420f4b53cSBarry Smith . sizes - array of length size (or `NULL`) providing the number of points in each partition 15520f4b53cSBarry Smith - points - array of length sum(sizes) (may be `NULL` iff sizes is `NULL`), a permutation of the points that groups those assigned to each partition in order (i.e., partition 0 first, partition 1 next, etc.) 156abe9303eSLisandro Dalcin 157abe9303eSLisandro Dalcin Level: developer 158abe9303eSLisandro Dalcin 15920f4b53cSBarry Smith Note: 160abe9303eSLisandro Dalcin It is safe to free the sizes and points arrays after use in this routine. 161abe9303eSLisandro Dalcin 162*a4e35b19SJacob Faibussowitsch .seealso: `DMPlexDistribute()`, `PetscPartitionerCreate()` 163abe9303eSLisandro Dalcin @*/ 164d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[]) 165d71ae5a4SJacob Faibussowitsch { 166abe9303eSLisandro Dalcin PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data; 167abe9303eSLisandro Dalcin PetscInt proc, numPoints; 168abe9303eSLisandro Dalcin 169abe9303eSLisandro Dalcin PetscFunctionBegin; 170abe9303eSLisandro Dalcin PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL); 1714f572ea9SToby Isaac if (sizes) PetscAssertPointer(sizes, 3); 1724f572ea9SToby Isaac if (points) PetscAssertPointer(points, 4); 1739566063dSJacob Faibussowitsch PetscCall(PetscSectionDestroy(&p->section)); 1749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&p->partition)); 1759566063dSJacob Faibussowitsch PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)part), &p->section)); 1769566063dSJacob Faibussowitsch PetscCall(PetscSectionSetChart(p->section, 0, size)); 177abe9303eSLisandro Dalcin if (sizes) { 17848a46eb9SPierre Jolivet for (proc = 0; proc < size; ++proc) PetscCall(PetscSectionSetDof(p->section, proc, sizes[proc])); 179abe9303eSLisandro Dalcin } 1809566063dSJacob Faibussowitsch PetscCall(PetscSectionSetUp(p->section)); 1819566063dSJacob Faibussowitsch PetscCall(PetscSectionGetStorageSize(p->section, &numPoints)); 1829566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)part), numPoints, points, PETSC_COPY_VALUES, &p->partition)); 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 184abe9303eSLisandro Dalcin } 185abe9303eSLisandro Dalcin 186abe9303eSLisandro Dalcin /*@ 187abe9303eSLisandro Dalcin PetscPartitionerShellSetRandom - Set the flag to use a random partition 188abe9303eSLisandro Dalcin 18920f4b53cSBarry Smith Collective 190abe9303eSLisandro Dalcin 191abe9303eSLisandro Dalcin Input Parameters: 19220f4b53cSBarry Smith + part - The `PetscPartitioner` 193abe9303eSLisandro Dalcin - random - The flag to use a random partition 194abe9303eSLisandro Dalcin 195abe9303eSLisandro Dalcin Level: intermediate 196abe9303eSLisandro Dalcin 197*a4e35b19SJacob Faibussowitsch .seealso: `PetscPartitionerShellGetRandom()`, `PetscPartitionerCreate()` 198abe9303eSLisandro Dalcin @*/ 199d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random) 200d71ae5a4SJacob Faibussowitsch { 201abe9303eSLisandro Dalcin PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data; 202abe9303eSLisandro Dalcin 203abe9303eSLisandro Dalcin PetscFunctionBegin; 204abe9303eSLisandro Dalcin PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL); 205abe9303eSLisandro Dalcin p->random = random; 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 207abe9303eSLisandro Dalcin } 208abe9303eSLisandro Dalcin 209abe9303eSLisandro Dalcin /*@ 210abe9303eSLisandro Dalcin PetscPartitionerShellGetRandom - get the flag to use a random partition 211abe9303eSLisandro Dalcin 21220f4b53cSBarry Smith Collective 213abe9303eSLisandro Dalcin 214abe9303eSLisandro Dalcin Input Parameter: 21520f4b53cSBarry Smith . part - The `PetscPartitioner` 216abe9303eSLisandro Dalcin 217abe9303eSLisandro Dalcin Output Parameter: 218abe9303eSLisandro Dalcin . random - The flag to use a random partition 219abe9303eSLisandro Dalcin 220abe9303eSLisandro Dalcin Level: intermediate 221abe9303eSLisandro Dalcin 222*a4e35b19SJacob Faibussowitsch .seealso: `PetscPartitionerShellSetRandom()`, `PetscPartitionerCreate()` 223abe9303eSLisandro Dalcin @*/ 224d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random) 225d71ae5a4SJacob Faibussowitsch { 226abe9303eSLisandro Dalcin PetscPartitioner_Shell *p = (PetscPartitioner_Shell *)part->data; 227abe9303eSLisandro Dalcin 228abe9303eSLisandro Dalcin PetscFunctionBegin; 229abe9303eSLisandro Dalcin PetscValidHeaderSpecificType(part, PETSCPARTITIONER_CLASSID, 1, PETSCPARTITIONERSHELL); 2304f572ea9SToby Isaac PetscAssertPointer(random, 2); 231abe9303eSLisandro Dalcin *random = p->random; 2323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 233abe9303eSLisandro Dalcin } 234