xref: /petsc/src/dm/partitioner/impls/gather/partgather.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
1abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h> /*I "petscpartitioner.h" I*/
2abe9303eSLisandro Dalcin 
3abe9303eSLisandro Dalcin typedef struct {
4abe9303eSLisandro Dalcin   PetscInt dummy;
5abe9303eSLisandro Dalcin } PetscPartitioner_Gather;
6abe9303eSLisandro Dalcin 
79371c9d4SSatish Balay static PetscErrorCode PetscPartitionerDestroy_Gather(PetscPartitioner part) {
8abe9303eSLisandro Dalcin   PetscFunctionBegin;
99566063dSJacob Faibussowitsch   PetscCall(PetscFree(part->data));
10abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
11abe9303eSLisandro Dalcin }
12abe9303eSLisandro Dalcin 
139371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_Gather_ASCII(PetscPartitioner part, PetscViewer viewer) {
14abe9303eSLisandro Dalcin   PetscFunctionBegin;
15abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
16abe9303eSLisandro Dalcin }
17abe9303eSLisandro Dalcin 
189371c9d4SSatish Balay static PetscErrorCode PetscPartitionerView_Gather(PetscPartitioner part, PetscViewer viewer) {
19abe9303eSLisandro Dalcin   PetscBool iascii;
20abe9303eSLisandro Dalcin 
21abe9303eSLisandro Dalcin   PetscFunctionBegin;
22abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
23abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
249566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
259566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscPartitionerView_Gather_ASCII(part, viewer));
26abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
27abe9303eSLisandro Dalcin }
28abe9303eSLisandro Dalcin 
299371c9d4SSatish Balay static PetscErrorCode PetscPartitionerPartition_Gather(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition) {
30abe9303eSLisandro Dalcin   PetscInt np;
31abe9303eSLisandro Dalcin 
32abe9303eSLisandro Dalcin   PetscFunctionBegin;
339566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition));
349566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetDof(partSection, 0, numVertices));
359566063dSJacob Faibussowitsch   for (np = 1; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, 0));
36abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
37abe9303eSLisandro Dalcin }
38abe9303eSLisandro Dalcin 
399371c9d4SSatish Balay static PetscErrorCode PetscPartitionerInitialize_Gather(PetscPartitioner part) {
40abe9303eSLisandro Dalcin   PetscFunctionBegin;
41abe9303eSLisandro Dalcin   part->noGraph        = PETSC_TRUE;
42abe9303eSLisandro Dalcin   part->ops->view      = PetscPartitionerView_Gather;
43abe9303eSLisandro Dalcin   part->ops->destroy   = PetscPartitionerDestroy_Gather;
44abe9303eSLisandro Dalcin   part->ops->partition = PetscPartitionerPartition_Gather;
45abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
46abe9303eSLisandro Dalcin }
47abe9303eSLisandro Dalcin 
48abe9303eSLisandro Dalcin /*MC
49abe9303eSLisandro Dalcin   PETSCPARTITIONERGATHER = "gather" - A PetscPartitioner object
50abe9303eSLisandro Dalcin 
51abe9303eSLisandro Dalcin   Level: intermediate
52abe9303eSLisandro Dalcin 
53db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
54abe9303eSLisandro Dalcin M*/
55abe9303eSLisandro Dalcin 
569371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Gather(PetscPartitioner part) {
57abe9303eSLisandro Dalcin   PetscPartitioner_Gather *p;
58abe9303eSLisandro Dalcin 
59abe9303eSLisandro Dalcin   PetscFunctionBegin;
60abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
61*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
62abe9303eSLisandro Dalcin   part->data = p;
63abe9303eSLisandro Dalcin 
649566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_Gather(part));
65abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
66abe9303eSLisandro Dalcin }
67