xref: /petsc/src/vec/is/sf/impls/basic/gather/sfgather.c (revision c73d2cf668e6709a5c9d6d1d49ca0c86e2ced9fc)
1 #include <../src/vec/is/sf/impls/basic/gatherv/sfgatherv.h>
2 
3 #define PetscSFPackGet_Gather PetscSFPackGet_Allgatherv
4 
5 /* Reuse the type. The difference is some fields (i.e., displs, recvcounts) are not used in Gather, which is not a big deal */
6 typedef PetscSF_Allgatherv PetscSF_Gather;
7 
8 PETSC_INTERN PetscErrorCode PetscSFBcastAndOpBegin_Gather(PetscSF sf,MPI_Datatype unit,PetscMemType rootmtype,const void *rootdata,PetscMemType leafmtype,void *leafdata,MPI_Op op)
9 {
10   PetscErrorCode       ierr;
11   PetscSFPack          link;
12   PetscMPIInt          rank,sendcount;
13   MPI_Comm             comm;
14   char                 *recvbuf;
15 
16   PetscFunctionBegin;
17   ierr = PetscSFPackGet_Gather(sf,unit,rootmtype,rootdata,leafmtype,leafdata,&link);CHKERRQ(ierr);
18   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
19 
20   if (op == MPIU_REPLACE) {
21     recvbuf = (char*)leafdata;
22   } else {
23     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
24     if (!link->leafbuf[leafmtype] && !rank) {ierr = PetscMallocWithMemType(leafmtype,sf->nleaves*link->unitbytes,(void**)&link->leafbuf[leafmtype]);CHKERRQ(ierr);}
25     recvbuf = link->leafbuf[leafmtype];
26   }
27 
28   ierr = PetscMPIIntCast(sf->nroots,&sendcount);CHKERRQ(ierr);
29   ierr = MPIU_Igather(rootdata,sendcount,unit,recvbuf,sendcount,unit,0/*rank 0*/,comm,link->rootreqs[PETSCSF_ROOT2LEAF_BCAST][rootmtype]);CHKERRQ(ierr);
30   PetscFunctionReturn(0);
31 }
32 
33 static PetscErrorCode PetscSFReduceBegin_Gather(PetscSF sf,MPI_Datatype unit,PetscMemType leafmtype,const void *leafdata,PetscMemType rootmtype,void *rootdata,MPI_Op op)
34 {
35   PetscErrorCode       ierr;
36   PetscSFPack          link;
37   PetscMPIInt          recvcount;
38   MPI_Comm             comm;
39   void                 *recvbuf;
40 
41   PetscFunctionBegin;
42   ierr = PetscSFPackGet_Gather(sf,unit,rootmtype,rootdata,leafmtype,leafdata,&link);CHKERRQ(ierr);
43   ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr);
44 
45   if (op == MPIU_REPLACE) {
46     recvbuf = (char*)rootdata;
47   } else {
48     if (!link->rootbuf[rootmtype]) {ierr = PetscMallocWithMemType(rootmtype,sf->nroots*link->unitbytes,(void**)&link->rootbuf[rootmtype]);CHKERRQ(ierr);}
49     recvbuf = link->rootbuf[rootmtype];
50   }
51 
52   ierr = PetscMPIIntCast(sf->nroots,&recvcount);CHKERRQ(ierr);
53   ierr = MPIU_Iscatter(leafdata,recvcount,unit,recvbuf,recvcount,unit,0/*rank 0*/,comm,link->rootreqs[PETSCSF_LEAF2ROOT_REDUCE][rootmtype]);CHKERRQ(ierr);
54   PetscFunctionReturn(0);
55 }
56 
57 PETSC_INTERN PetscErrorCode PetscSFCreate_Gather(PetscSF sf)
58 {
59   PetscErrorCode  ierr;
60   PetscSF_Gather  *dat = (PetscSF_Gather*)sf->data;
61 
62   PetscFunctionBegin;
63   /* Inherit from Allgatherv */
64   sf->ops->Reset           = PetscSFReset_Allgatherv;
65   sf->ops->Destroy         = PetscSFDestroy_Allgatherv;
66   sf->ops->GetGraph        = PetscSFGetGraph_Allgatherv;
67   sf->ops->GetRootRanks    = PetscSFGetRootRanks_Allgatherv;
68   sf->ops->GetLeafRanks    = PetscSFGetLeafRanks_Allgatherv;
69   sf->ops->BcastAndOpEnd   = PetscSFBcastAndOpEnd_Allgatherv;
70   sf->ops->ReduceEnd       = PetscSFReduceEnd_Allgatherv;
71   sf->ops->FetchAndOpEnd   = PetscSFFetchAndOpEnd_Allgatherv;
72   sf->ops->CreateLocalSF   = PetscSFCreateLocalSF_Allgatherv;
73 
74   /* Inherit from Gatherv */
75   sf->ops->FetchAndOpBegin = PetscSFFetchAndOpBegin_Gatherv;
76 
77   /* Gather stuff */
78   sf->ops->BcastAndOpBegin = PetscSFBcastAndOpBegin_Gather;
79   sf->ops->ReduceBegin     = PetscSFReduceBegin_Gather;
80 
81   ierr     = PetscNewLog(sf,&dat);CHKERRQ(ierr);
82   sf->data = (void*)dat;
83   PetscFunctionReturn(0);
84 }
85 
86