From: Anthony Geay Date: Tue, 10 Mar 2020 20:08:22 +0000 (+0100) Subject: Preparation for integration X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=3f0e7ce788a7306a2855a9d6653fe53d9597c8ef;p=tools%2Fmedcoupling.git Preparation for integration --- diff --git a/src/ParaMEDMEM/CommInterface.hxx b/src/ParaMEDMEM/CommInterface.hxx index 861013366..32ff28dde 100644 --- a/src/ParaMEDMEM/CommInterface.hxx +++ b/src/ParaMEDMEM/CommInterface.hxx @@ -79,10 +79,12 @@ namespace MEDCoupling int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) const { return MPI_Allgather(sendbuf,sendcount, sendtype, recvbuf, recvcount, recvtype, comm); } + int allGatherV(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int recvcounts[], + const int displs[], MPI_Datatype recvtype, MPI_Comm comm) { return MPI_Allgatherv(sendbuf,sendcount,sendtype,recvbuf,recvcounts,displs,recvtype,comm); } int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype, void* recvbuf, int recvcount, MPI_Datatype recvtype, MPI_Comm comm) const { return MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm); } - int allToAllV(void* sendbuf, int* sendcounts, int* senddispls, + int allToAllV(const void* sendbuf, int* sendcounts, int* senddispls, MPI_Datatype sendtype, void* recvbuf, int* recvcounts, int* recvdispls, MPI_Datatype recvtype, MPI_Comm comm) const { return MPI_Alltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm); } @@ -91,7 +93,8 @@ namespace MEDCoupling MPI_Op op, int root, MPI_Comm comm) const { return MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm); } int allReduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) const { return MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm); } public: - static std::unique_ptr ToIntArray(const std::unique_ptr& arr, std::size_t size) + template + static std::unique_ptr ToIntArray(const std::unique_ptr& arr, std::size_t size) { std::unique_ptr ret(new int[size]); std::copy(arr.get(),arr.get()+size,ret.get()); diff --git a/src/ParaMEDMEM/ParaUMesh.cxx b/src/ParaMEDMEM/ParaUMesh.cxx index c03a6cd23..1d4970f00 100644 --- a/src/ParaMEDMEM/ParaUMesh.cxx +++ b/src/ParaMEDMEM/ParaUMesh.cxx @@ -53,23 +53,28 @@ ParaUMesh::ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, Dat throw INTERP_KERNEL::Exception("ParaUMesh constructor : mismatch between # cells and len of global # cells."); } +/*! +* This method computes the cells part of distributed mesh lying on \a globalNodeIds nodes. +* The input \a globalNodeIds are not supposed to reside on the current process. +*/ MCAuto ParaUMesh::getCellIdsLyingOnNodes(DataArrayIdType *globalNodeIds, bool fullyIn) const { if(fullyIn) throw INTERP_KERNEL::Exception("ParaUMesh::getCellIdsLyingOnNodes : not implemented yet for fullyIn == True !"); MPI_Comm comm(MPI_COMM_WORLD); + CommInterface ci; int size; - MPI_Comm_size(comm,&size); - std::unique_ptr nbOfElems(new long[size]),nbOfElems2(new long[size]); - long nbOfNodeIdsLoc(FromIdType(globalNodeIds->getNumberOfTuples())); - MPI_Allgather(&nbOfNodeIdsLoc,1,MPI_LONG,nbOfElems.get(),1,MPI_LONG,comm); - long nbOfNodeIdsSum(std::accumulate(nbOfElems.get(),nbOfElems.get()+size,0L)); + ci.commSize(comm,&size); + std::unique_ptr nbOfElems(new mcIdType[size]),nbOfElems2(new mcIdType[size]); + mcIdType nbOfNodeIdsLoc(globalNodeIds->getNumberOfTuples()); + ci.allGather(&nbOfNodeIdsLoc,1,MPI_ID_TYPE,nbOfElems.get(),1,MPI_ID_TYPE,comm); + mcIdType nbOfNodeIdsSum(std::accumulate(nbOfElems.get(),nbOfElems.get()+size,0)); std::vector< MCAuto > tabs(size); { - std::unique_ptr allGlobalNodeIds(new int[nbOfNodeIdsSum]); - std::unique_ptr nbOfElemsInt( CommInterface::ToIntArray(nbOfElems,size) ); + std::unique_ptr allGlobalNodeIds(new mcIdType[nbOfNodeIdsSum]); + std::unique_ptr nbOfElemsInt( CommInterface::ToIntArray(nbOfElems,size) ); std::unique_ptr offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ); - MPI_Allgatherv(globalNodeIds->begin(),globalNodeIds->getNumberOfTuples(),MPI_INT,allGlobalNodeIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_INT,comm); + ci.allGatherV(globalNodeIds->begin(),globalNodeIds->getNumberOfTuples(),MPI_ID_TYPE,allGlobalNodeIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm); long offset(0); for(int curRk = 0 ; curRk < size ; ++curRk) { @@ -80,21 +85,21 @@ MCAuto ParaUMesh::getCellIdsLyingOnNodes(DataArrayIdType *globa MCAuto localNodeIdsToLocate(_node_global->findIdForEach(globalNodeIdsCaptured->begin(),globalNodeIdsCaptured->end())); MCAuto localCellCaptured(_mesh->getCellIdsLyingOnNodes(localNodeIdsToLocate->begin(),localNodeIdsToLocate->end(),false)); MCAuto localCellCapturedGlob(_cell_global->selectByTupleIdSafe(localCellCaptured->begin(),localCellCaptured->end())); - nbOfElems[curRk] = FromIdType(localCellCapturedGlob->getNumberOfTuples()); + nbOfElems[curRk] = localCellCapturedGlob->getNumberOfTuples(); tabs[curRk] = localCellCapturedGlob; } } std::vector tabss(tabs.begin(),tabs.end()); MCAuto cells(DataArrayIdType::Aggregate(tabss)); - MPI_Alltoall(nbOfElems.get(),1,MPI_LONG,nbOfElems2.get(),1,MPI_LONG,comm); + ci.allToAll(nbOfElems.get(),1,MPI_ID_TYPE,nbOfElems2.get(),1,MPI_ID_TYPE,comm); long nbOfCellIdsSum(std::accumulate(nbOfElems2.get(),nbOfElems2.get()+size,0L)); MCAuto cellIdsFromProcs(DataArrayIdType::New()); cellIdsFromProcs->alloc(nbOfCellIdsSum,1); { - std::unique_ptr nbOfElemsInt( CommInterface::ToIntArray(nbOfElems,size) ),nbOfElemsOutInt( CommInterface::ToIntArray(nbOfElems2,size) ); + std::unique_ptr nbOfElemsInt( CommInterface::ToIntArray(nbOfElems,size) ),nbOfElemsOutInt( CommInterface::ToIntArray(nbOfElems2,size) ); std::unique_ptr offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ), offsetsOut( CommInterface::ComputeOffset(nbOfElemsOutInt,size) ); - MPI_Alltoallv(cells->begin(),nbOfElemsInt.get(),offsetsIn.get(),MPI_INT, - cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_INT,comm); + ci.allToAllV(cells->begin(),nbOfElemsInt.get(),offsetsIn.get(),MPI_INT, + cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_INT,comm); } cellIdsFromProcs->sort(); return cellIdsFromProcs;