xref: /petsc/src/dm/partitioner/impls/gather/partgather.c (revision abe9303e6325b68e0d8957680d58b261e7a295d5)
1*abe9303eSLisandro Dalcin #include <petsc/private/partitionerimpl.h>        /*I "petscpartitioner.h" I*/
2*abe9303eSLisandro Dalcin 
3*abe9303eSLisandro Dalcin typedef struct {
4*abe9303eSLisandro Dalcin   PetscInt dummy;
5*abe9303eSLisandro Dalcin } PetscPartitioner_Gather;
6*abe9303eSLisandro Dalcin 
7*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerDestroy_Gather(PetscPartitioner part)
8*abe9303eSLisandro Dalcin {
9*abe9303eSLisandro Dalcin   PetscErrorCode ierr;
10*abe9303eSLisandro Dalcin 
11*abe9303eSLisandro Dalcin   PetscFunctionBegin;
12*abe9303eSLisandro Dalcin   ierr = PetscFree(part->data);CHKERRQ(ierr);
13*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
14*abe9303eSLisandro Dalcin }
15*abe9303eSLisandro Dalcin 
16*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Gather_ASCII(PetscPartitioner part, PetscViewer viewer)
17*abe9303eSLisandro Dalcin {
18*abe9303eSLisandro Dalcin   PetscFunctionBegin;
19*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
20*abe9303eSLisandro Dalcin }
21*abe9303eSLisandro Dalcin 
22*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerView_Gather(PetscPartitioner part, PetscViewer viewer)
23*abe9303eSLisandro Dalcin {
24*abe9303eSLisandro Dalcin   PetscBool      iascii;
25*abe9303eSLisandro Dalcin   PetscErrorCode ierr;
26*abe9303eSLisandro Dalcin 
27*abe9303eSLisandro Dalcin   PetscFunctionBegin;
28*abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
29*abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
30*abe9303eSLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
31*abe9303eSLisandro Dalcin   if (iascii) {ierr = PetscPartitionerView_Gather_ASCII(part, viewer);CHKERRQ(ierr);}
32*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
33*abe9303eSLisandro Dalcin }
34*abe9303eSLisandro Dalcin 
35*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerPartition_Gather(PetscPartitioner part, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection vertSection, PetscSection targetSection, PetscSection partSection, IS *partition)
36*abe9303eSLisandro Dalcin {
37*abe9303eSLisandro Dalcin   PetscInt       np;
38*abe9303eSLisandro Dalcin   PetscErrorCode ierr;
39*abe9303eSLisandro Dalcin 
40*abe9303eSLisandro Dalcin   PetscFunctionBegin;
41*abe9303eSLisandro Dalcin   ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr);
42*abe9303eSLisandro Dalcin   ierr = PetscSectionSetDof(partSection,0,numVertices);CHKERRQ(ierr);
43*abe9303eSLisandro Dalcin   for (np = 1; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, 0);CHKERRQ(ierr);}
44*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
45*abe9303eSLisandro Dalcin }
46*abe9303eSLisandro Dalcin 
47*abe9303eSLisandro Dalcin static PetscErrorCode PetscPartitionerInitialize_Gather(PetscPartitioner part)
48*abe9303eSLisandro Dalcin {
49*abe9303eSLisandro Dalcin   PetscFunctionBegin;
50*abe9303eSLisandro Dalcin   part->noGraph        = PETSC_TRUE;
51*abe9303eSLisandro Dalcin   part->ops->view      = PetscPartitionerView_Gather;
52*abe9303eSLisandro Dalcin   part->ops->destroy   = PetscPartitionerDestroy_Gather;
53*abe9303eSLisandro Dalcin   part->ops->partition = PetscPartitionerPartition_Gather;
54*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
55*abe9303eSLisandro Dalcin }
56*abe9303eSLisandro Dalcin 
57*abe9303eSLisandro Dalcin /*MC
58*abe9303eSLisandro Dalcin   PETSCPARTITIONERGATHER = "gather" - A PetscPartitioner object
59*abe9303eSLisandro Dalcin 
60*abe9303eSLisandro Dalcin   Level: intermediate
61*abe9303eSLisandro Dalcin 
62*abe9303eSLisandro Dalcin .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
63*abe9303eSLisandro Dalcin M*/
64*abe9303eSLisandro Dalcin 
65*abe9303eSLisandro Dalcin PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Gather(PetscPartitioner part)
66*abe9303eSLisandro Dalcin {
67*abe9303eSLisandro Dalcin   PetscPartitioner_Gather *p;
68*abe9303eSLisandro Dalcin   PetscErrorCode           ierr;
69*abe9303eSLisandro Dalcin 
70*abe9303eSLisandro Dalcin   PetscFunctionBegin;
71*abe9303eSLisandro Dalcin   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
72*abe9303eSLisandro Dalcin   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
73*abe9303eSLisandro Dalcin   part->data = p;
74*abe9303eSLisandro Dalcin 
75*abe9303eSLisandro Dalcin   ierr = PetscPartitionerInitialize_Gather(part);CHKERRQ(ierr);
76*abe9303eSLisandro Dalcin   PetscFunctionReturn(0);
77*abe9303eSLisandro Dalcin }
78*abe9303eSLisandro Dalcin 
79*abe9303eSLisandro Dalcin 
80