/*------------------------------- phasicFlow --------------------------------- O C enter of O O E ngineering and O O M ultiscale modeling of OOOOOOO F luid flow ------------------------------------------------------------------------------ Copyright (C): www.cemf.ir email: hamid.r.norouzi AT gmail.com ------------------------------------------------------------------------------ Licence: This file is part of phasicFlow code. It is a free software for simulating granular and multiphase flows. You can redistribute it and/or modify it under the terms of GNU General Public License v3 or any other later versions. phasicFlow is distributed to help others in their research in the field of granular and multiphase flows, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. -----------------------------------------------------------------------------*/ #ifndef __mpiCommunication_H__ #define __mpiCommunication_H__ #include "mpiTypes.hpp" #include "types.hpp" #include "span.hpp" #ifdef pFlow_Build_MPI namespace pFlow::MPI { extern DataType realx3Type__; extern DataType realx4Type__; extern DataType int32x3Type__; template auto constexpr Type() { return MPI_BYTE; } template auto constexpr sFactor() { return sizeof(T); } template auto constexpr Type() { return MPI_CHAR; } template auto constexpr sFactor() { return 1; } template auto constexpr Type() { return MPI_SHORT; } template auto constexpr sFactor() { return 1; } template auto constexpr Type() { return MPI_UNSIGNED_SHORT; } template auto constexpr sFactor() { return 1; } template auto constexpr Type() { return MPI_INT; } template auto constexpr sFactor() { return 1; } template<> auto constexpr Type() { return MPI_UNSIGNED; } template<> auto constexpr sFactor() { return 1; } template<> auto constexpr Type() { return MPI_LONG; } template<> auto constexpr sFactor() { return 1; } template<> auto constexpr Type() { return MPI_UNSIGNED_LONG; } template<> auto constexpr sFactor() { return 1; } template<> auto constexpr Type() { return MPI_FLOAT; } template<> auto constexpr sFactor() { return 1; } template<> auto constexpr Type() { return MPI_DOUBLE; } template<> auto constexpr sFactor() { return 1; } template<> inline auto Type() { return realx3Type__; } template<> auto constexpr sFactor() { return 1; } template<> inline auto Type() { return realx4Type__; } template<> auto constexpr sFactor() { return 1; } template<> inline auto Type() { return int32x3Type__; } template<> auto constexpr sFactor() { return 1; } /*inline auto createByteSequence(int sizeOfElement) { DataType newType; MPI_Type_contiguous(sizeOfElement, MPI_CHAR, &newType); MPI_Type_commit(&newType); return newType; }*/ inline auto TypeCommit(DataType* type) { return MPI_Type_commit(type); } inline auto TypeFree(DataType* type) { return MPI_Type_free(type); } template inline auto getCount(Status* status, int& count) { int lCount; auto res = MPI_Get_count(status, Type(), &lCount); count = lCount/sFactor(); return res; } template inline int convertIndex(const int& ind) { return ind*sFactor(); } template inline auto send(span data, int dest, int tag, Comm comm) { return MPI_Send( data.data(), sFactor()*data().size(), Type(), dest, tag, comm); } template inline auto recv(span data, int source, int tag, Comm comm, Status *status) { return MPI_Recv( data.data(), sFactor()*data.size(), Type(), source, tag, comm, status); } template inline auto scan(T sData, T& rData, Comm comm, Operation op = SumOp) { return MPI_Scan(&sData, &rData, sFactor()*1, Type(), op , comm ); } // gathering one scalar data to root processor template inline auto gather(T sendData, span& recvData, int root, Comm comm) { return MPI_Gather( &sendData, sFactor()*1, Type(), recvData.data(), sFactor()*1, Type(), root, comm); } template inline auto allGather(T sendData, span& recvData, Comm comm) { return MPI_Allgather( &sendData, sFactor()*1, Type(), recvData.data(), sFactor()*1, Type(), comm); } template inline auto scatter(span sendData, T& recvData, int root, Comm comm) { return MPI_Scatter( sendData.data(), sFactor()*1, Type(), &recvData, sFactor()*1, Type(), root, comm); } template inline auto Bcast(T& sendData, int root, Comm comm) { return MPI_Bcast( &sendData, sFactor()*1, Type(), root, comm); } template bool typeCreateIndexedBlock( span index, DataType &newType) { auto res = MPI_Type_create_indexed_block( index.size(), sFactor(), index.data(), Type(), &newType); if(res == Success) { TypeCommit(&newType); } else { return false; } return true; } template inline auto Gatherv ( span sendData, span& recvData, span recvCounts, span displs, int root, Comm comm) { return MPI_Gatherv( sendData.data(), sendData.size()*sFactor(), Type(), recvData.data(), recvCounts.data(), displs.data(), Type(), root, comm ); } inline auto Wait(Request* request, Status* status) { return MPI_Wait(request, status); } inline auto typeFree(DataType& type) { return MPI_Type_free(&type); } } #endif //pFlow_Build_MPI #endif //__mpiCommunication_H__