openPMD-api
ChunkInfo.hpp
1 /* Copyright 2020-2025 Franz Poeschel, Axel Huebl
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 #pragma once
22 
23 #include "openPMD/config.hpp"
24 
25 #include "openPMD/Dataset.hpp" // Offset, Extent
26 #include "openPMD/auxiliary/BlockSlicer.hpp"
27 
28 #if openPMD_HAVE_MPI
29 #include <mpi.h>
30 #endif
31 
32 #include <map>
33 #include <memory>
34 #include <string>
35 #include <vector>
36 
37 namespace openPMD
38 {
44 struct ChunkInfo
45 {
46  Offset offset;
47  Extent extent;
48 
49  /*
50  * If rank is smaller than zero, will be converted to zero.
51  */
52  explicit ChunkInfo() = default;
53  ChunkInfo(Offset, Extent);
54 
55  bool operator==(ChunkInfo const &other) const;
56 };
57 
72 {
73  unsigned int sourceID = 0;
74 
75  explicit WrittenChunkInfo() = default;
76  /*
77  * If rank is smaller than zero, will be converted to zero.
78  */
79  WrittenChunkInfo(Offset, Extent, int sourceID);
80  WrittenChunkInfo(Offset, Extent);
81 
82  bool operator==(WrittenChunkInfo const &other) const;
83 };
84 
85 using ChunkTable = std::vector<WrittenChunkInfo>;
86 
87 namespace chunk_assignment
88 {
92  using RankMeta = std::map<unsigned int, std::string>;
93 
104  using Assignment = std::map<unsigned int, std::vector<WrittenChunkInfo>>;
105 
115  template <typename Chunk_t>
116  void mergeChunks(std::vector<Chunk_t> &chunks);
117 
125  auto
126  mergeChunksFromSameSourceID(std::vector<WrittenChunkInfo> const &chunks)
127  -> std::map<unsigned int, std::vector<ChunkInfo>>;
128 
137  {
138  ChunkTable notAssigned;
139  Assignment assigned;
140 
141  explicit PartialAssignment() = default;
142  PartialAssignment(ChunkTable notAssigned);
143  PartialAssignment(ChunkTable notAssigned, Assignment assigned);
144  };
145 
154  struct Strategy
155  {
180  Assignment assign(
181  ChunkTable chunkTable,
182  RankMeta const &in,
183  RankMeta const &out,
184  size_t my_rank,
185  size_t num_ranks);
213  virtual Assignment assign(
214  PartialAssignment partialAssignment,
215  RankMeta const &in,
216  RankMeta const &out,
217  size_t my_rank,
218  size_t num_ranks) = 0;
219 
220  virtual std::unique_ptr<Strategy> clone() const = 0;
221 
222  virtual ~Strategy() = default;
223  };
224 
239  {
267  ChunkTable table,
268  RankMeta const &in,
269  RankMeta const &out,
270  size_t my_rank,
271  size_t num_ranks);
301  PartialAssignment partialAssignment,
302  RankMeta const &in,
303  RankMeta const &out,
304  size_t my_rank,
305  size_t num_ranks) = 0;
306 
307  virtual std::unique_ptr<PartialStrategy> clone() const = 0;
308 
309  virtual ~PartialStrategy() = default;
310  };
311 
325  {
327  std::unique_ptr<PartialStrategy> firstPass,
328  std::unique_ptr<Strategy> secondPass);
329 
330  virtual Assignment assign(
332  RankMeta const &in,
333  RankMeta const &out,
334  size_t my_rank,
335  size_t num_ranks) override;
336 
337  virtual std::unique_ptr<Strategy> clone() const override;
338 
339  private:
340  std::unique_ptr<PartialStrategy> m_firstPass;
341  std::unique_ptr<Strategy> m_secondPass;
342  };
343 
350  {
351  Assignment assign(
353  RankMeta const &in,
354  RankMeta const &out,
355  size_t my_rank,
356  size_t num_ranks) override;
357 
358  virtual std::unique_ptr<Strategy> clone() const override;
359  };
360 
368  {
369  Assignment assign(
371  RankMeta const &in,
372  RankMeta const &out,
373  size_t my_rank,
374  size_t num_ranks) override;
375 
376  virtual std::unique_ptr<Strategy> clone() const override;
377  };
378 
389  struct Blocks : Strategy
390  {
391  Assignment assign(
393  RankMeta const &in,
394  RankMeta const &out,
395  size_t my_rank,
396  size_t num_ranks) override;
397 
398  [[nodiscard]] std::unique_ptr<Strategy> clone() const override;
399  };
400 
408  {
409  Assignment assign(
411  RankMeta const &in,
412  RankMeta const &out,
413  size_t my_rank,
414  size_t num_ranks) override;
415 
416  [[nodiscard]] std::unique_ptr<Strategy> clone() const override;
417  };
418 
428  {
429  ByHostname(std::unique_ptr<Strategy> withinNode);
430 
433  RankMeta const &in,
434  RankMeta const &out,
435  size_t my_rank,
436  size_t num_ranks) override;
437 
438  virtual std::unique_ptr<PartialStrategy> clone() const override;
439 
440  private:
441  std::unique_ptr<Strategy> m_withinNode;
442  };
443 
455  {
457  std::unique_ptr<auxiliary::BlockSlicer> blockSlicer,
458  Extent totalExtent);
459 
460  Assignment assign(
462  RankMeta const &in,
463  RankMeta const &out,
464  size_t my_rank,
465  size_t num_ranks) override;
466 
467  virtual std::unique_ptr<Strategy> clone() const override;
468 
469  private:
470  std::unique_ptr<auxiliary::BlockSlicer> blockSlicer;
471  Extent totalExtent;
472  };
473 
487  {
488  size_t splitAlongDimension = 0;
489 
494  BinPacking(size_t splitAlongDimension = 0);
495 
496  Assignment assign(
498  RankMeta const &in,
499  RankMeta const &out,
500  size_t my_rank,
501  size_t num_ranks) override;
502 
503  virtual std::unique_ptr<Strategy> clone() const override;
504  };
505 
516  {
517  explicit FailingStrategy();
518 
519  Assignment assign(
521  RankMeta const &in,
522  RankMeta const &out,
523  size_t my_rank,
524  size_t num_ranks) override;
525 
526  virtual std::unique_ptr<Strategy> clone() const override;
527  };
528 
539  {
540  explicit DiscardingStrategy();
541 
542  Assignment assign(
544  RankMeta const &in,
545  RankMeta const &out,
546  size_t my_rank,
547  size_t num_ranks) override;
548 
549  virtual std::unique_ptr<Strategy> clone() const override;
550  };
551 } // namespace chunk_assignment
552 
553 namespace host_info
554 {
560  enum class Method
561  {
562  POSIX_HOSTNAME,
563  MPI_PROCESSOR_NAME
564  };
565 
572  bool methodAvailable(Method);
573 
581  std::string byMethod(Method);
582 
583 #if openPMD_HAVE_MPI
594  chunk_assignment::RankMeta byMethodCollective(MPI_Comm, Method);
595 #endif
596 } // namespace host_info
597 } // namespace openPMD
598 
599 #undef openPMD_POSIX_AVAILABLE
Public definitions of openPMD-api.
Definition: Date.cpp:29
Represents the meta info around a chunk in a dataset.
Definition: ChunkInfo.hpp:45
Extent extent
size of the chunk
Definition: ChunkInfo.hpp:47
Offset offset
origin of the chunk
Definition: ChunkInfo.hpp:46
Represents the meta info around a chunk that has been written by some data producing application.
Definition: ChunkInfo.hpp:72
unsigned int sourceID
ID of the data source containing the chunk.
Definition: ChunkInfo.hpp:73
Strategy that tries to assign chunks in a balanced manner without arbitrarily cutting chunks.
Definition: ChunkInfo.hpp:487
BinPacking(size_t splitAlongDimension=0)
Definition: ChunkInfo.cpp:736
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:740
Alternative to RoundRobin, but instead gives every reader a sequential range of blocks.
Definition: ChunkInfo.hpp:390
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:366
Blocks at processs level.
Definition: ChunkInfo.hpp:408
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:390
Slice the n-dimensional dataset into hyperslabs and distribute chunks according to them.
Definition: ChunkInfo.hpp:455
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:701
Strategy that assigns chunks to be read by processes within the same host that produced the chunk.
Definition: ChunkInfo.hpp:428
PartialAssignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:434
Strategy that purposefully discards leftover chunk from the PartialAssignment.
Definition: ChunkInfo.hpp:539
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:872
Strategy that purposefully fails when the PartialAssignment has leftover chunks.
Definition: ChunkInfo.hpp:516
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:847
Combine a PartialStrategy and a Strategy to obtain a Strategy working in two phases.
Definition: ChunkInfo.hpp:325
virtual Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:274
Return type for partial chunk assignment strategies.
Definition: ChunkInfo.hpp:137
A chunk distribution strategy that guarantees no complete distribution.
Definition: ChunkInfo.hpp:239
PartialAssignment assign(ChunkTable table, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks)
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:253
virtual PartialAssignment assign(PartialAssignment partialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks)=0
Assign chunks to be loaded to reading processes.
Simple strategy that assigns produced chunks to reading processes in a round-Robin manner.
Definition: ChunkInfo.hpp:350
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:296
Round-Robin at process level.
Definition: ChunkInfo.hpp:368
Assignment assign(PartialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks) override
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:329
Interface for a chunk distribution strategy.
Definition: ChunkInfo.hpp:155
virtual Assignment assign(PartialAssignment partialAssignment, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks)=0
Assign chunks to be loaded to reading processes.
Assignment assign(ChunkTable chunkTable, RankMeta const &in, RankMeta const &out, size_t my_rank, size_t num_ranks)
Assign chunks to be loaded to reading processes.
Definition: ChunkInfo.cpp:224