xref: /petsc/src/dm/partitioner/impls/gather/partgather.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 
7d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerDestroy_Gather(PetscPartitioner part)
8d71ae5a4SJacob Faibussowitsch {
9abe9303eSLisandro Dalcin   PetscFunctionBegin;
109566063dSJacob Faibussowitsch   PetscCall(PetscFree(part->data));
11*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12abe9303eSLisandro Dalcin }
13abe9303eSLisandro Dalcin 
14d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Gather_ASCII(PetscPartitioner part, PetscViewer viewer)
15d71ae5a4SJacob Faibussowitsch {
16abe9303eSLisandro Dalcin   PetscFunctionBegin;
17*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
18abe9303eSLisandro Dalcin }
19abe9303eSLisandro Dalcin 
20d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerView_Gather(PetscPartitioner part, PetscViewer viewer)
21d71ae5a4SJacob Faibussowitsch {
22abe9303eSLisandro Dalcin   PetscBool iascii;
23abe9303eSLisandro Dalcin 
24abe9303eSLisandro Dalcin   PetscFunctionBegin;
25abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
26abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
289566063dSJacob Faibussowitsch   if (iascii) PetscCall(PetscPartitionerView_Gather_ASCII(part, viewer));
29*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30abe9303eSLisandro Dalcin }
31abe9303eSLisandro Dalcin 
32d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerPartition_Gather(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
33d71ae5a4SJacob Faibussowitsch {
34abe9303eSLisandro Dalcin   PetscInt np;
35abe9303eSLisandro Dalcin 
36abe9303eSLisandro Dalcin   PetscFunctionBegin;
379566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition));
389566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetDof(partSection, 0, numVertices));
399566063dSJacob Faibussowitsch   for (np = 1; np < nparts; ++np) PetscCall(PetscSectionSetDof(partSection, np, 0));
40*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41abe9303eSLisandro Dalcin }
42abe9303eSLisandro Dalcin 
43d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscPartitionerInitialize_Gather(PetscPartitioner part)
44d71ae5a4SJacob Faibussowitsch {
45abe9303eSLisandro Dalcin   PetscFunctionBegin;
46abe9303eSLisandro Dalcin   part->noGraph        = PETSC_TRUE;
47abe9303eSLisandro Dalcin   part->ops->view      = PetscPartitionerView_Gather;
48abe9303eSLisandro Dalcin   part->ops->destroy   = PetscPartitionerDestroy_Gather;
49abe9303eSLisandro Dalcin   part->ops->partition = PetscPartitionerPartition_Gather;
50*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51abe9303eSLisandro Dalcin }
52abe9303eSLisandro Dalcin 
53abe9303eSLisandro Dalcin /*MC
54abe9303eSLisandro Dalcin   PETSCPARTITIONERGATHER = "gather" - A PetscPartitioner object
55abe9303eSLisandro Dalcin 
56abe9303eSLisandro Dalcin   Level: intermediate
57abe9303eSLisandro Dalcin 
58db781477SPatrick Sanan .seealso: `PetscPartitionerType`, `PetscPartitionerCreate()`, `PetscPartitionerSetType()`
59abe9303eSLisandro Dalcin M*/
60abe9303eSLisandro Dalcin 
61d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Gather(PetscPartitioner part)
62d71ae5a4SJacob Faibussowitsch {
63abe9303eSLisandro Dalcin   PetscPartitioner_Gather *p;
64abe9303eSLisandro Dalcin 
65abe9303eSLisandro Dalcin   PetscFunctionBegin;
66abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
674dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&p));
68abe9303eSLisandro Dalcin   part->data = p;
69abe9303eSLisandro Dalcin 
709566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerInitialize_Gather(part));
71*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
72abe9303eSLisandro Dalcin }
73