KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
Classes | Namespaces | Macros
data_communicator.h File Reference
#include <string>
#include <iostream>
#include <type_traits>
#include "containers/array_1d.h"
#include "containers/flags.h"
#include "includes/define.h"
#include "includes/mpi_serializer.h"
Include dependency graph for data_communicator.h:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Classes

class  Kratos::DataCommunicator
 Serial (do-nothing) version of a wrapper class for MPI communication. More...
 

Namespaces

 Kratos
 REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
 

Macros

#define KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(Size1, Size2, CheckedFunction)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SYNC_SHAPE_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_LOC_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE(...)
 Wrappers for MPI_Scatter and MPI_Scatterv calls. More...
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE(...)
 Wrappers for MPI_Gather, MPI_Gatherv and MPI_Allgather calls. More...
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_PUBLIC_INTERFACE_FOR_TYPE(...)
 
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_IMPLEMENTATION_FOR_TYPE(...)
 

Functions

Input and output
std::istream & Kratos::operator>> (std::istream &rIStream, DataCommunicator &rThis)
 input stream function More...
 
std::ostream & Kratos::operator<< (std::ostream &rOStream, const DataCommunicator &rThis)
 output stream function More...
 

Macro Definition Documentation

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE (   ...)
Value:
virtual __VA_ARGS__ SumAll(const __VA_ARGS__& rLocalValue) const { return rLocalValue; } \
virtual std::vector<__VA_ARGS__> SumAll(const std::vector<__VA_ARGS__>& rLocalValues) const { \
return rLocalValues; \
} \
virtual void SumAll(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rLocalValues.size(), rGlobalValues.size(), "SumAll"); \
rGlobalValues = SumAll(rLocalValues); \
} \
virtual __VA_ARGS__ MinAll(const __VA_ARGS__& rLocalValue) const { return rLocalValue; } \
virtual std::vector<__VA_ARGS__> MinAll(const std::vector<__VA_ARGS__>& rLocalValues) const { \
return rLocalValues; \
} \
virtual void MinAll(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rLocalValues.size(), rGlobalValues.size(), "MinAll"); \
rGlobalValues = MinAll(rLocalValues); \
} \
virtual __VA_ARGS__ MaxAll(const __VA_ARGS__& rLocalValue) const { return rLocalValue; } \
virtual std::vector<__VA_ARGS__> MaxAll(const std::vector<__VA_ARGS__>& rLocalValues) const { \
return rLocalValues; \
} \
virtual void MaxAll(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rLocalValues.size(), rGlobalValues.size(), "MaxAll"); \
rGlobalValues = MaxAll(rLocalValues); \
}

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_LOC_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_LOC_INTERFACE_FOR_TYPE (   ...)
Value:
virtual std::pair<__VA_ARGS__, int> MinLocAll(const __VA_ARGS__& rLocalValue) const { return std::pair<__VA_ARGS__, int>(rLocalValue, 0); } \
virtual std::pair<__VA_ARGS__, int> MaxLocAll(const __VA_ARGS__& rLocalValue) const { return std::pair<__VA_ARGS__, int>(rLocalValue, 0); }

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE (   ...)
Value:
virtual void BroadcastImpl(__VA_ARGS__& rBuffer, const int SourceRank) const {} \
virtual void BroadcastImpl(std::vector<__VA_ARGS__>& rBuffer, const int SourceRank) const {} \

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE (   ...)

Wrappers for MPI_Gather, MPI_Gatherv and MPI_Allgather calls.

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_IMPLEMENTATION_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_IMPLEMENTATION_FOR_TYPE (   ...)
Value:
KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE(__VA_ARGS__) \
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE(...)
Definition: data_communicator.h:141

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_PUBLIC_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_PUBLIC_INTERFACE_FOR_TYPE (   ...)
Value:
KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE(__VA_ARGS__) \
KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE(__VA_ARGS__) \
KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE(__VA_ARGS__) \
KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE(__VA_ARGS__) \
KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SYNC_SHAPE_INTERFACE_FOR_TYPE(__VA_ARGS__) \
#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE(...)
Definition: data_communicator.h:50

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE (   ...)
Value:
virtual __VA_ARGS__ Sum(const __VA_ARGS__& rLocalValue, const int Root) const { return rLocalValue; } \
virtual std::vector<__VA_ARGS__> Sum(const std::vector<__VA_ARGS__>& rLocalValues, const int Root) const { \
return rLocalValues; \
} \
virtual void Sum(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues, const int Root) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rLocalValues.size(), rGlobalValues.size(), "Sum"); \
rGlobalValues = Sum(rLocalValues, Root); \
} \
virtual __VA_ARGS__ Min(const __VA_ARGS__& rLocalValue, const int Root) const { return rLocalValue; } \
virtual std::vector<__VA_ARGS__> Min(const std::vector<__VA_ARGS__>& rLocalValues, const int Root) const { \
return rLocalValues; \
} \
virtual void Min(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues, const int Root) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rLocalValues.size(), rGlobalValues.size(), "Min"); \
rGlobalValues = Min(rLocalValues, Root); \
} \
virtual __VA_ARGS__ Max(const __VA_ARGS__& rLocalValue, const int Root) const { return rLocalValue; } \
virtual std::vector<__VA_ARGS__> Max(const std::vector<__VA_ARGS__>& rLocalValues, const int Root) const { \
return rLocalValues; \
} \
virtual void Max(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues, const int Root) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rLocalValues.size(), rGlobalValues.size(), "Max"); \
rGlobalValues = Max(rLocalValues, Root); \
} \

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE (   ...)
Value:
virtual __VA_ARGS__ ScanSum(const __VA_ARGS__& rLocalValue) const { return rLocalValue; } \
virtual std::vector<__VA_ARGS__> ScanSum(const std::vector<__VA_ARGS__>& rLocalValues) const { \
return rLocalValues; \
} \
virtual void ScanSum(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rPartialSums) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rLocalValues.size(), rPartialSums.size(), "ScanSum"); \
rPartialSums = ScanSum(rLocalValues); \
} \

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE (   ...)
Value:
virtual std::vector<__VA_ARGS__> Scatter(const std::vector<__VA_ARGS__>& rSendValues, const int SourceRank) const { \
KRATOS_ERROR_IF( Rank() != SourceRank ) \
<< "Communication between different ranks is not possible with a serial DataCommunicator." << std::endl; \
return rSendValues; \
} \
virtual void Scatter( \
const std::vector<__VA_ARGS__>& rSendValues, std::vector<__VA_ARGS__>& rRecvValues, const int SourceRank) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rSendValues.size(),rRecvValues.size(),"Scatter"); \
rRecvValues = Scatter(rSendValues, SourceRank); \
} \
virtual std::vector<__VA_ARGS__> Scatterv(const std::vector<std::vector<__VA_ARGS__>>& rSendValues, const int SourceRank) const { \
KRATOS_ERROR_IF( Rank() != SourceRank ) \
<< "Communication between different ranks is not possible with a serial DataCommunicator." << std::endl; \
KRATOS_ERROR_IF( static_cast<unsigned int>(Size()) != rSendValues.size() ) \
<< "Unexpected number of sends in DataCommuncatior::Scatterv (serial DataCommunicator always assumes a single process)." << std::endl; \
return rSendValues[0]; \
} \
virtual void Scatterv( \
const std::vector<__VA_ARGS__>& rSendValues, const std::vector<int>& rSendCounts, const std::vector<int>& rSendOffsets, \
std::vector<__VA_ARGS__>& rRecvValues, const int SourceRank) const { \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rRecvValues.size(), rSendValues.size(), "Scatterv (values check)"); \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rSendCounts.size(), 1, "Scatterv (counts check)"); \
KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK(rSendOffsets.size(), 1, "Scatterv (offsets check)"); \
KRATOS_ERROR_IF( Rank() != SourceRank ) \
<< "Communication between different ranks is not possible with a serial DataCommunicator." << std::endl; \
rRecvValues = rSendValues; \
} \
TSpaceType::IndexType Size(TSpaceType &dummy, typename TSpaceType::VectorType const &rV)
Definition: add_strategies_to_python.cpp:111

Wrappers for MPI_Scatter and MPI_Scatterv calls.

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE (   ...)

◆ KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SYNC_SHAPE_INTERFACE_FOR_TYPE

#define KRATOS_BASE_DATA_COMMUNICATOR_DECLARE_SYNC_SHAPE_INTERFACE_FOR_TYPE (   ...)
Value:
virtual bool SynchronizeShape(__VA_ARGS__&) const { return false; } \
virtual bool SynchronizeShape( \
const __VA_ARGS__& rSendValue, const int SendDestination, const int SendTag, \
__VA_ARGS__& rRecvValue, const int RecvSource, const int RecvTag) const { return false; } \

◆ KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK

#define KRATOS_DATA_COMMUNICATOR_DEBUG_SIZE_CHECK (   Size1,
  Size2,
  CheckedFunction 
)
Value:
<< "Input error in call to DataCommunicator::" << CheckedFunction \
<< ": The sizes of the local and distributed buffers do not match." << std::endl;
#define KRATOS_DEBUG_ERROR_IF(conditional)
Definition: exception.h:171
TSpaceType::IndexType Size1(TSpaceType &dummy, typename TSpaceType::MatrixType const &rM)
Definition: add_strategies_to_python.cpp:117
TSpaceType::IndexType Size2(TSpaceType &dummy, typename TSpaceType::MatrixType const &rM)
Definition: add_strategies_to_python.cpp:123