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