24 #include "openPMD/config.hpp" 27 #include "openPMD/Datatype.hpp" 28 #include "openPMD/Series.hpp" 43 template<
typename Duration >
46 MPI_Comm communicator;
61 typename decltype( Series::iterations )::key_type
94 std::string compression,
96 std::string extension,
99 typename decltype( Series::iterations )::key_type iterations,
122 std::string compression,
124 std::string extension,
127 typename decltype( Series::iterations)::key_type iterations
140 template<
typename Dummy >
146 MPI_Datatype dt = MPI_CHAR;
148 template<
typename Dummy >
154 MPI_Datatype dt = MPI_UNSIGNED_CHAR;
156 template<
typename Dummy >
162 MPI_Datatype dt = MPI_SHORT;
164 template<
typename Dummy >
170 MPI_Datatype dt = MPI_INT;
172 template<
typename Dummy >
178 MPI_Datatype dt = MPI_LONG;
180 template<
typename Dummy >
186 MPI_Datatype dt = MPI_FLOAT;
188 template<
typename Dummy >
194 MPI_Datatype dt = MPI_DOUBLE;
196 template<
typename Dummy >
202 MPI_Datatype dt = MPI_UNSIGNED_SHORT;
204 template<
typename Dummy >
210 MPI_Datatype dt = MPI_UNSIGNED;
212 template<
typename Dummy >
218 MPI_Datatype dt = MPI_UNSIGNED_LONG;
220 template<
typename Dummy >
226 MPI_Datatype dt = MPI_LONG_DOUBLE;
228 template<
typename Dummy >
234 MPI_Datatype dt = MPI_LONG_LONG_INT;
237 MPIDatatype< typename Duration::rep > m_mpiDatatype;
238 MPI_Datatype mpiType = m_mpiDatatype.dt;
244 template<
typename Duration >
247 std::string compression,
249 std::string extension,
252 typename decltype( Series::iterations )::key_type iterations,
259 using rep =
typename Duration::rep;
274 rank < threadSize ? 0 : MPI_UNDEFINED,
279 if( rank < threadSize )
288 rep * recv =
nullptr;
289 if( rank == rootThread )
291 recv =
new rep[2 * threadSize];
294 if( restricted != MPI_COMM_NULL )
309 if( rank == rootThread )
311 for(
int i = 0; i < threadSize; i++ )
313 Duration dWrite { recv[2 * i] };
314 Duration dRead { recv[2 * i + 1] };
334 if( restricted != MPI_COMM_NULL )
336 MPI_Comm_free( &restricted );
340 template<
typename Duration >
345 template<
typename Duration >
351 std::string compression,
353 std::string extension,
356 typename decltype( Series::iterations )::key_type iterations
377 throw std::runtime_error(
"Requested report not found. (Reports are available on the root thread only)" );
void addReport(int rootThread, std::string compression, uint8_t level, std::string extension, int threadSize, Datatype dt, typename decltype(Series::iterations)::key_type iterations, std::pair< Duration, Duration > const &report)
Add results for a certain compression strategy and level.
Definition: MPIBenchmarkReport.hpp:245
std::map< std::tuple< int, std::string, uint8_t, std::string, int, Datatype, typename decltype(Series::iterations)::key_type >, std::pair< Duration, Duration > > durations
Time needed for writing and reading per compression strategy and level.
Definition: MPIBenchmarkReport.hpp:67
Datatype
Concrete datatype of an object available at runtime.
Definition: Datatype.hpp:42
The report for a single benchmark produced by <openPMD/benchmark/mpi/MPIBenchmark>.
Definition: MPIBenchmarkReport.hpp:44
Public definitions of openPMD-api.
Definition: Date.cpp:29
std::pair< Duration, Duration > getReport(int rank, std::string compression, uint8_t level, std::string extension, int threadSize, Datatype dt, typename decltype(Series::iterations)::key_type iterations)
Retrieve the time measured for a certain compression strategy.
Definition: MPIBenchmarkReport.hpp:349