openPMD-api
MPIBenchmark.hpp
1 /* Copyright 2018-2020 Franz Poeschel
2  *
3  * This file is part of openPMD-api.
4  *
5  * openPMD-api is free software: you can redistribute it and/or modify
6  * it under the terms of of either the GNU General Public License or
7  * the GNU Lesser General Public License as published by
8  * the Free Software Foundation, either version 3 of the License, or
9  * (at your option) any later version.
10  *
11  * openPMD-api is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License and the GNU Lesser General Public License
15  * for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * and the GNU Lesser General Public License along with openPMD-api.
19  * If not, see <http://www.gnu.org/licenses/>.
20  */
21 
22 #pragma once
23 
24 #include "openPMD/config.hpp"
25 #if openPMD_HAVE_MPI
26 
27 #include "RandomDatasetFiller.hpp"
28 
29 #include "openPMD/openPMD.hpp"
30 #include "openPMD/benchmark/mpi/MPIBenchmarkReport.hpp"
31 #include "openPMD/benchmark/mpi/DatasetFiller.hpp"
32 #include "openPMD/benchmark/mpi/BlockSlicer.hpp"
33 
34 #include <mpi.h>
35 
36 #include <chrono>
37 #include <exception>
38 #include <iostream>
39 #include <sstream>
40 #include <utility>
41 #include <vector>
42 #include <tuple>
43 #include <set>
44 
45 
46 namespace openPMD
47 {
55  template< typename DatasetFillerProvider >
57  {
58 
59  public:
60  using extentT = Extent::value_type;
61  MPI_Comm communicator = MPI_COMM_WORLD;
62 
66  Extent totalExtent;
67 
68  std::shared_ptr< BlockSlicer > m_blockSlicer;
69 
70  DatasetFillerProvider m_dfp;
71 
72 
88  std::string basePath,
89  Extent tExtent,
90  std::shared_ptr< BlockSlicer > blockSlicer,
91  DatasetFillerProvider dfp,
92  MPI_Comm comm = MPI_COMM_WORLD
93  );
94 
105  void addConfiguration(
106  std::string compression,
107  uint8_t compressionLevel,
108  std::string backend,
109  Datatype dt,
110  typename decltype( Series::iterations )::key_type iterations,
111  int threadSize
112  );
113 
125  void addConfiguration(
126  std::string compression,
127  uint8_t compressionLevel,
128  std::string backend,
129  Datatype dt,
130  typename decltype( Series::iterations)::key_type iterations
131  );
132 
133  void resetConfigurations( );
134 
135 
144  template< typename Clock >
146  int rootThread = 0
147  );
148 
149  private:
150  std::string m_basePath;
151  std::vector<
152  std::tuple<
153  std::string,
154  uint8_t,
155  std::string,
156  int,
157  Datatype,
158  typename decltype( Series::iterations)::key_type>>
159  m_configurations;
160 
161  enum Config
162  {
163  COMPRESSION = 0,
164  COMPRESSION_LEVEL,
165  BACKEND,
166  NRANKS,
167  DTYPE,
168  ITERATIONS
169  };
170 
171  std::pair<
172  Offset,
173  Extent
174  > slice( int size );
175 
182  template< typename Clock >
183  struct BenchmarkExecution
184  {
186 
187 
188  explicit BenchmarkExecution( MPIBenchmark< DatasetFillerProvider > * benchmark ) :
189  m_benchmark { benchmark }
190  {}
191 
192 
205  template<
206  typename T
207  >
208  typename Clock::duration writeBenchmark(
209  std::string const & compression,
210  uint8_t level,
211  Offset & offset,
212  Extent & extent,
213  std::string const & extension,
214  std::shared_ptr< DatasetFiller< T >> datasetFiller,
215  typename decltype( Series::iterations)::key_type iterations
216  );
217 
227  template<
228  typename T
229  >
230  typename Clock::duration readBenchmark(
231  Offset & offset,
232  Extent & extent,
233  std::string extension,
234  typename decltype( Series::iterations)::key_type iterations
235  );
236 
237  template< typename T >
238  void operator()(
240  int rootThread = 0
241  );
242 
243 
244  template< int n >
245  void operator()(
247  int
248  );
249  };
250  };
251 
252 
253  // Implementation
254 
255 
256 
257 
258 
259  template< typename DatasetFillerProvider >
260  template< typename Clock >
263  int rootThread
264  )
265  {
266  MPIBenchmarkReport< typename Clock::duration > res{this->communicator};
267  BenchmarkExecution< Clock > exec { this };
268 
269  std::set< Datatype > datatypes;
270  for( auto const & conf: m_configurations )
271  {
272  datatypes.insert( std::get< DTYPE >( conf ) );
273  }
274  for( Datatype dt: datatypes )
275  {
276  switchType(
277  dt,
278  exec,
279  res,
280  rootThread
281  );
282  }
283 
284  return res;
285  }
286 
287 
288  template< typename DatasetFillerProvider >
290  std::string basePath,
291  Extent tExtent,
292  std::shared_ptr< BlockSlicer > blockSlicer,
293  DatasetFillerProvider dfp,
294  MPI_Comm comm
295  ):
296  communicator { comm },
297  totalExtent { std::move( tExtent ) },
298  m_blockSlicer { std::move( blockSlicer ) },
299  m_dfp { dfp },
300  m_basePath { std::move( basePath ) }
301  {
302  if( m_blockSlicer == nullptr )
303  throw std::runtime_error("Argument blockSlicer cannot be a nullptr!");
304  }
305 
306 
307  template< typename DatasetFillerProvider >
308  std::pair<
309  Offset,
310  Extent
312  {
313  int actualSize;
314  MPI_Comm_size(
315  this->communicator,
316  &actualSize
317  );
318  int rank;
319  MPI_Comm_rank(
320  this->communicator,
321  &rank
322  );
323  size = std::min(
324  size,
325  actualSize
326  );
327  return m_blockSlicer->sliceBlock(
328  totalExtent,
329  size,
330  rank
331  );
332  }
333 
334 
335  template< typename DatasetFillerProvider >
337  std::string compression,
338  uint8_t compressionLevel,
339  std::string backend,
340  Datatype dt,
341  typename decltype( Series::iterations)::key_type iterations,
342  int threadSize
343  )
344  {
345  this->m_configurations
346  .emplace_back(
347  compression,
348  compressionLevel,
349  backend,
350  threadSize,
351  dt,
352  iterations
353  );
354  }
355 
356 
357  template< typename DatasetFillerProvider >
359  std::string compression,
360  uint8_t compressionLevel,
361  std::string backend,
362  Datatype dt,
363  typename decltype( Series::iterations)::key_type iterations
364  )
365  {
366  int size;
367  MPI_Comm_size(
368  communicator,
369  &size
370  );
372  compression,
373  compressionLevel,
374  backend,
375  dt,
376  iterations,
377  size
378  );
379  }
380 
381 
382  template< typename DatasetFillerProvider >
384  {
385  this->m_compressions
386  .clear( );
387  }
388 
389 
390  template< typename DatasetFillerProvider >
391  template< typename Clock >
392  template< typename T >
393  typename Clock::duration
395  std::string const & compression,
396  uint8_t level,
397  Offset & offset,
398  Extent & extent,
399  std::string const & extension,
400  std::shared_ptr< DatasetFiller< T >> datasetFiller,
401  typename decltype( Series::iterations)::key_type iterations
402  )
403  {
404  MPI_Barrier( m_benchmark->communicator );
405  auto start = Clock::now( );
406 
407  // open file for writing
408  Series series = Series(
409  m_benchmark->m_basePath + "." + extension,
411  m_benchmark->communicator
412  );
413 
414  for( typename decltype( Series::iterations)::key_type i = 0;
415  i < iterations;
416  i++ )
417  {
418  auto writeData = datasetFiller->produceData( );
419 
420 
422  id =
423  series.iterations[i].meshes["id"][MeshRecordComponent::SCALAR];
424 
425  Datatype datatype = determineDatatype( writeData );
426  Dataset dataset = Dataset(
427  datatype,
428  m_benchmark->totalExtent
429  );
430  if( !compression.empty( ) )
431  {
432  dataset.setCompression(
433  compression,
434  level
435  );
436  }
437 
438  id.resetDataset( dataset );
439 
440  series.flush( );
441 
442  id.storeChunk< T >(
443  writeData,
444  offset,
445  extent
446  );
447  series.flush( );
448  }
449 
450  MPI_Barrier( m_benchmark->communicator );
451  auto end = Clock::now( );
452 
453  // deduct the time needed for data generation
454  for( typename decltype( Series::iterations)::key_type i = 0;
455  i < iterations;
456  i++ )
457  {
458  datasetFiller->produceData( );
459  }
460  auto deduct = Clock::now( );
461 
462  return end - start - ( deduct - end );
463  }
464 
465 
466  template< typename DatasetFillerProvider >
467  template< typename Clock >
468  template< typename T >
469  typename Clock::duration
471  Offset & offset,
472  Extent & extent,
473  std::string extension,
474  typename decltype( Series::iterations)::key_type iterations
475  )
476  {
477  MPI_Barrier( m_benchmark->communicator );
478  // let every thread measure time
479  auto start = Clock::now( );
480 
481  Series series = Series(
482  m_benchmark->m_basePath + "." + extension,
484  m_benchmark->communicator
485  );
486 
487  for( typename decltype( Series::iterations)::key_type i = 0;
488  i < iterations;
489  i++ )
490  {
492  id =
493  series.iterations[i].meshes["id"][MeshRecordComponent::SCALAR];
494 
495 
496  auto chunk_data = id.loadChunk< T >(
497  offset,
498  extent
499  );
500  series.flush( );
501  }
502 
503  MPI_Barrier( m_benchmark->communicator );
504  auto end = Clock::now( );
505  return end - start;
506  }
507 
508 
509  template< typename DatasetFillerProvider >
510  template< typename Clock >
511  template< typename T >
512  void
515  int rootThread
516  )
517  {
518  Datatype dt = determineDatatype< T >( );
519  auto dsf = std::dynamic_pointer_cast< DatasetFiller< T>>(
520  m_benchmark->m_dfp
521  .template operator(
522  )< T >( )
523  );
524  for( auto const & config: m_benchmark->m_configurations )
525  {
526  std::string compression;
527  uint8_t compressionLevel;
528  std::string backend;
529  int size;
530  Datatype dt2;
531  typename decltype( Series::iterations)::key_type iterations;
532  std::tie(
533  compression,
534  compressionLevel,
535  backend,
536  size,
537  dt2,
538  iterations
539  ) = config;
540 
541  if( dt != dt2 )
542  {
543  continue;
544  }
545 
546  auto localCuboid = m_benchmark->slice( size );
547 
548  extentT blockSize = 1;
549  for( auto ext: localCuboid.second )
550  {
551  blockSize *= ext;
552  }
553  dsf->setNumberOfItems( blockSize );
554 
555  auto writeTime = writeBenchmark< T >(
556  compression,
557  compressionLevel,
558  localCuboid.first,
559  localCuboid.second,
560  backend,
561  dsf,
562  iterations
563  );
564  auto readTime = readBenchmark< T >(
565  localCuboid.first,
566  localCuboid.second,
567  backend,
568  iterations
569  );
570  report.addReport(
571  rootThread,
572  compression,
573  compressionLevel,
574  backend,
575  size,
576  dt2,
577  iterations,
578  std::make_pair(
579  writeTime,
580  readTime
581  )
582  );
583 
584  }
585  }
586 
587 
588  template< typename DatasetFillerProvider >
589  template< typename Clock >
590  template< int n >
591  void
594  int
595  )
596  {
597  throw std::runtime_error( "Unknown/unsupported datatype requested to be benchmarked." );
598  }
599 
600 }
601 
602 #endif
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
Definition: Dataset.hpp:36
An abstract class to create one iteration of data per thread.
Definition: DatasetFiller.hpp:37
Extent totalExtent
Total extent of the hypercuboid used in the benchmark.
Definition: MPIBenchmark.hpp:66
create new series and truncate existing (files)
MPIBenchmarkReport< typename Clock::duration > runBenchmark(int rootThread=0)
Main function for running a benchmark.
Definition: MPIBenchmark.hpp:262
Datatype
Concrete datatype of an object available at runtime.
Definition: Datatype.hpp:42
Root level of the openPMD hierarchy.
Definition: Series.hpp:64
Class representing a benchmark.
Definition: MPIBenchmark.hpp:56
The report for a single benchmark produced by <openPMD/benchmark/mpi/MPIBenchmark>.
Definition: MPIBenchmarkReport.hpp:44
void flush()
Execute all required remaining IO operations to write or read data.
Definition: Series.cpp:359
MPIBenchmark(std::string basePath, Extent tExtent, std::shared_ptr< BlockSlicer > blockSlicer, DatasetFillerProvider dfp, MPI_Comm comm=MPI_COMM_WORLD)
Construct an MPI benchmark manually.
Definition: MPIBenchmark.hpp:289
ReturnType switchType(Datatype dt, Action action, Args &&... args)
Generalizes switching over an openPMD datatype.
Definition: Datatype.hpp:642
Public definitions of openPMD-api.
Definition: Date.cpp:29
Definition: MeshRecordComponent.hpp:30
open series as read-only, fails if series is not found
void addConfiguration(std::string compression, uint8_t compressionLevel, std::string backend, Datatype dt, typename decltype(Series::iterations)::key_type iterations, int threadSize)
Definition: MPIBenchmark.hpp:336