Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

enable msvc for dataman and table #2422

Merged
merged 6 commits into from
Aug 12, 2020
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
24 changes: 10 additions & 14 deletions cmake/DetectOptions.cmake
Original file line number Diff line number Diff line change
Expand Up @@ -166,13 +166,11 @@ if(ZeroMQ_FOUND)
endif()

# DataMan
if(NOT MSVC)
if(ZeroMQ_FOUND)
if(ADIOS2_USE_DataMan STREQUAL AUTO)
set(ADIOS2_HAVE_DataMan TRUE)
elseif(ADIOS2_USE_DataMan)
set(ADIOS2_HAVE_DataMan TRUE)
endif()
if(ZeroMQ_FOUND)
if(ADIOS2_USE_DataMan STREQUAL AUTO)
set(ADIOS2_HAVE_DataMan TRUE)
elseif(ADIOS2_USE_DataMan)
set(ADIOS2_HAVE_DataMan TRUE)
endif()
endif()

Expand All @@ -188,13 +186,11 @@ if(NOT MSVC)
endif()

# Table
if(NOT MSVC)
if(ZeroMQ_FOUND)
if(ADIOS2_USE_Table STREQUAL AUTO)
set(ADIOS2_HAVE_Table TRUE)
elseif(ADIOS2_USE_Table)
set(ADIOS2_HAVE_Table TRUE)
endif()
if(ZeroMQ_FOUND)
if(ADIOS2_USE_Table STREQUAL AUTO)
set(ADIOS2_HAVE_Table TRUE)
elseif(ADIOS2_USE_Table)
set(ADIOS2_HAVE_Table TRUE)
endif()
endif()

Expand Down
9 changes: 5 additions & 4 deletions source/adios2/engine/ssc/SscHelper.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -348,7 +348,8 @@ void MPI_Gatherv64(const void *sendbuf, uint64_t sendcount,
}
}

MPI_Waitall(requests.size(), requests.data(), MPI_STATUSES_IGNORE);
MPI_Waitall(static_cast<int>(requests.size()), requests.data(),
MPI_STATUSES_IGNORE);
}

void MPI_Gatherv64OneSidedPull(const void *sendbuf, uint64_t sendcount,
Expand Down Expand Up @@ -447,9 +448,9 @@ void MPI_Gatherv64OneSidedPush(const void *sendbuf, uint64_t sendcount,
{
MPI_Put(reinterpret_cast<const char *>(sendbuf) +
(sendcount - sendcountvar) * sendTypeSize,
sendcountvar, sendtype, root,
displs[mpiRank] + sendcount - sendcountvar, sendcountvar,
sendtype, win);
static_cast<int>(sendcountvar), sendtype, root,
displs[mpiRank] + sendcount - sendcountvar,
static_cast<int>(sendcountvar), sendtype, win);
sendcountvar = 0;
}
}
Expand Down
33 changes: 18 additions & 15 deletions source/adios2/engine/ssc/SscReader.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -13,6 +13,7 @@
#include "adios2/helper/adiosCommMPI.h"
#include "adios2/helper/adiosFunctions.h"
#include "adios2/helper/adiosJSONcomplex.h"
#include "adios2/helper/adiosMpiHandshake.h"
#include "nlohmann/json.hpp"

namespace adios2
Expand Down Expand Up @@ -77,8 +78,8 @@ StepStatus SscReader::BeginStep(const StepMode stepMode,
{
if (m_MpiMode == "twosided")
{
MPI_Status statuses[m_MpiRequests.size()];
MPI_Waitall(m_MpiRequests.size(), m_MpiRequests.data(), statuses);
MPI_Waitall(m_MpiRequests.size(), m_MpiRequests.data(),
MPI_STATUS_IGNORE);
m_MpiRequests.clear();
}
else if (m_MpiMode == "onesidedfencepush")
Expand Down Expand Up @@ -190,9 +191,9 @@ void SscReader::EndStep()
for (const auto &i : m_AllReceivingWriterRanks)
{
m_MpiRequests.emplace_back();
MPI_Irecv(m_Buffer.data() + i.second.first, i.second.second,
MPI_CHAR, i.first, 0, m_StreamComm,
&m_MpiRequests.back());
MPI_Irecv(m_Buffer.data() + i.second.first,
static_cast<int>(i.second.second), MPI_CHAR, i.first,
0, m_StreamComm, &m_MpiRequests.back());
}
}
else if (m_MpiMode == "onesidedfencepush")
Expand All @@ -208,19 +209,19 @@ void SscReader::EndStep()
MPI_Win_fence(0, m_MpiWin);
for (const auto &i : m_AllReceivingWriterRanks)
{
MPI_Get(m_Buffer.data() + i.second.first, i.second.second,
MPI_CHAR, i.first, 0, i.second.second, MPI_CHAR,
m_MpiWin);
MPI_Get(m_Buffer.data() + i.second.first,
static_cast<int>(i.second.second), MPI_CHAR, i.first, 0,
static_cast<int>(i.second.second), MPI_CHAR, m_MpiWin);
}
}
else if (m_MpiMode == "onesidedpostpull")
{
MPI_Win_start(m_MpiAllWritersGroup, 0, m_MpiWin);
for (const auto &i : m_AllReceivingWriterRanks)
{
MPI_Get(m_Buffer.data() + i.second.first, i.second.second,
MPI_CHAR, i.first, 0, i.second.second, MPI_CHAR,
m_MpiWin);
MPI_Get(m_Buffer.data() + i.second.first,
static_cast<int>(i.second.second), MPI_CHAR, i.first, 0,
static_cast<int>(i.second.second), MPI_CHAR, m_MpiWin);
}
}
}
Expand Down Expand Up @@ -267,8 +268,9 @@ bool SscReader::SyncWritePattern()
m_StreamComm);
std::vector<char> localVec(maxLocalSize, '\0');
std::vector<char> globalVec(maxLocalSize * m_StreamSize);
MPI_Allgather(localVec.data(), maxLocalSize, MPI_CHAR, globalVec.data(),
maxLocalSize, MPI_CHAR, m_StreamComm);
MPI_Allgather(localVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
globalVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
m_StreamComm);

ssc::LocalJsonToGlobalJson(globalVec, maxLocalSize, m_StreamSize,
m_GlobalWritePatternJson);
Expand Down Expand Up @@ -416,8 +418,9 @@ void SscReader::SyncReadPattern()
std::vector<char> localVec(maxLocalSize, '\0');
std::memcpy(localVec.data(), localStr.c_str(), localStr.size());
std::vector<char> globalVec(maxLocalSize * m_StreamSize);
MPI_Allgather(localVec.data(), maxLocalSize, MPI_CHAR, globalVec.data(),
maxLocalSize, MPI_CHAR, m_StreamComm);
MPI_Allgather(localVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
globalVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
m_StreamComm);

// deserialize global metadata Json
nlohmann::json globalJson;
Expand Down
9 changes: 4 additions & 5 deletions source/adios2/engine/ssc/SscReader.tcc
Original file line number Diff line number Diff line change
Expand Up @@ -67,17 +67,16 @@ void SscReader::GetDeferredCommon(Variable<std::string> &variable,
for (const auto &i : m_AllReceivingWriterRanks)
{
MPI_Win_lock(MPI_LOCK_SHARED, i.first, 0, m_MpiWin);
MPI_Get(m_Buffer.data() + i.second.first, i.second.second,
MPI_CHAR, i.first, 0, i.second.second, MPI_CHAR,
m_MpiWin);
MPI_Get(m_Buffer.data() + i.second.first,
static_cast<int>(i.second.second), MPI_CHAR, i.first, 0,
static_cast<int>(i.second.second), MPI_CHAR, m_MpiWin);
MPI_Win_unlock(i.first, m_MpiWin);
}
}

for (const auto &i : m_AllReceivingWriterRanks)
{
const auto &v = m_GlobalWritePattern[i.first];
for (const auto &b : v)
for (const auto &b : m_GlobalWritePattern[i.first])
{
if (b.name == variable.m_Name)
{
Expand Down
32 changes: 18 additions & 14 deletions source/adios2/engine/ssc/SscWriter.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -109,7 +109,8 @@ void SscWriter::EndStep()
for (const auto &i : m_AllSendingReaderRanks)
{
m_MpiRequests.emplace_back();
MPI_Isend(m_Buffer.data(), m_Buffer.size(), MPI_CHAR,
MPI_Isend(m_Buffer.data(),
static_cast<int>(m_Buffer.size()), MPI_CHAR,
i.first, 0, m_StreamComm, &m_MpiRequests.back());
}
}
Expand All @@ -118,8 +119,9 @@ void SscWriter::EndStep()
MPI_Win_fence(0, m_MpiWin);
for (const auto &i : m_AllSendingReaderRanks)
{
MPI_Put(m_Buffer.data(), m_Buffer.size(), MPI_CHAR, i.first,
i.second.first, m_Buffer.size(), MPI_CHAR,
MPI_Put(m_Buffer.data(), static_cast<int>(m_Buffer.size()),
MPI_CHAR, i.first, i.second.first,
static_cast<int>(m_Buffer.size()), MPI_CHAR,
m_MpiWin);
}
}
Expand All @@ -128,8 +130,9 @@ void SscWriter::EndStep()
MPI_Win_start(m_MpiAllReadersGroup, 0, m_MpiWin);
for (const auto &i : m_AllSendingReaderRanks)
{
MPI_Put(m_Buffer.data(), m_Buffer.size(), MPI_CHAR, i.first,
i.second.first, m_Buffer.size(), MPI_CHAR,
MPI_Put(m_Buffer.data(), static_cast<int>(m_Buffer.size()),
MPI_CHAR, i.first, i.second.first,
static_cast<int>(m_Buffer.size()), MPI_CHAR,
m_MpiWin);
}
}
Expand Down Expand Up @@ -160,8 +163,8 @@ void SscWriter::MpiWait()
{
if (m_MpiMode == "twosided")
{
MPI_Status statuses[m_MpiRequests.size()];
MPI_Waitall(m_MpiRequests.size(), m_MpiRequests.data(), statuses);
MPI_Waitall(static_cast<int>(m_MpiRequests.size()),
m_MpiRequests.data(), MPI_STATUSES_IGNORE);
m_MpiRequests.clear();
}
else if (m_MpiMode == "onesidedfencepush")
Expand Down Expand Up @@ -230,8 +233,9 @@ void SscWriter::SyncWritePattern(bool finalStep)
std::vector<char> localVec(maxLocalSize, '\0');
std::memcpy(localVec.data(), localStr.data(), localStr.size());
std::vector<char> globalVec(maxLocalSize * m_StreamSize, '\0');
MPI_Allgather(localVec.data(), maxLocalSize, MPI_CHAR, globalVec.data(),
maxLocalSize, MPI_CHAR, m_StreamComm);
MPI_Allgather(localVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
globalVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
m_StreamComm);

// deserialize global metadata Json
nlohmann::json globalJson;
Expand All @@ -258,8 +262,9 @@ void SscWriter::SyncReadPattern()
m_StreamComm);
std::vector<char> localVec(maxLocalSize, '\0');
std::vector<char> globalVec(maxLocalSize * m_StreamSize);
MPI_Allgather(localVec.data(), maxLocalSize, MPI_CHAR, globalVec.data(),
maxLocalSize, MPI_CHAR, m_StreamComm);
MPI_Allgather(localVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
globalVec.data(), static_cast<int>(maxLocalSize), MPI_CHAR,
m_StreamComm);

// deserialize global metadata Json
nlohmann::json globalJson;
Expand Down Expand Up @@ -321,7 +326,7 @@ void SscWriter::CalculatePosition(ssc::BlockVecVec &writerVecVec,
auto currentReaderOverlapWriterRanks =
CalculateOverlap(writerVecVec, readerRankMap);
size_t bufferPosition = 0;
for (size_t rank = 0; rank < writerVecVec.size(); ++rank)
for (int rank = 0; rank < writerVecVec.size(); ++rank)
{
bool hasOverlap = false;
for (const auto r : currentReaderOverlapWriterRanks)
Expand Down Expand Up @@ -382,8 +387,7 @@ void SscWriter::DoClose(const int transportIndex)
MPI_Isend(m_Buffer.data(), 1, MPI_CHAR, i.first, 0,
m_StreamComm, &requests.back());
}
MPI_Status statuses[requests.size()];
MPI_Waitall(requests.size(), requests.data(), statuses);
MPI_Waitall(requests.size(), requests.data(), MPI_STATUS_IGNORE);
}
else if (m_MpiMode == "onesidedfencepush")
{
Expand Down