Serial Examples
The serial API provides sequential, one-process read and write access. Most users will use this for exploration and processing of their data.
Reading
C++
#include <openPMD/openPMD.hpp>
#include <cstddef>
#include <iostream>
#include <memory>
using std::cout;
using namespace openPMD;
int main()
{
Series series =
Series("../samples/git-sample/data%T.h5", Access::READ_ONLY);
cout << "Read a Series with openPMD standard version " << series.openPMD()
<< '\n';
cout << "The Series contains " << series.snapshots().size()
<< " iterations:";
for (auto const &i : series.snapshots())
cout << "\n\t" << i.first;
cout << '\n';
Iteration i = series.snapshots()[100];
cout << "Iteration 100 contains " << i.meshes.size() << " meshes:";
for (auto const &m : i.meshes)
cout << "\n\t" << m.first;
cout << '\n';
cout << "Iteration 100 contains " << i.particles.size()
<< " particle species:";
for (auto const &ps : i.particles)
{
cout << "\n\t" << ps.first;
for (auto const &r : ps.second)
{
cout << "\n\t" << r.first;
cout << '\n';
}
}
openPMD::ParticleSpecies electrons = i.particles["electrons"];
std::shared_ptr<double> charge = electrons["charge"].loadChunk<double>();
series.flush();
cout << "And the first electron particle has a charge = "
<< charge.get()[0];
cout << '\n';
MeshRecordComponent E_x = i.meshes["E"]["x"];
Extent extent = E_x.getExtent();
cout << "Field E/x has shape (";
for (auto const &dim : extent)
cout << dim << ',';
cout << ") and has datatype " << E_x.getDatatype() << '\n';
Offset chunk_offset = {1, 1, 1};
Extent chunk_extent = {2, 2, 1};
// Loading without explicit datatype here
auto chunk_data = E_x.loadChunkVariant(chunk_offset, chunk_extent);
cout << "Queued the loading of a single chunk from disk, "
"ready to execute\n";
series.flush();
cout << "Chunk has been read from disk\n"
<< "Read chunk contains:\n";
std::visit(
[&chunk_offset, &chunk_extent](auto &shared_ptr) {
for (size_t row = 0; row < chunk_extent[0]; ++row)
{
for (size_t col = 0; col < chunk_extent[1]; ++col)
cout << "\t" << '(' << row + chunk_offset[0] << '|'
<< col + chunk_offset[1] << '|' << 1 << ")\t"
<< shared_ptr.get()[row * chunk_extent[1] + col];
cout << '\n';
}
},
chunk_data);
auto all_data = E_x.loadChunk<double>();
// The iteration can be closed in order to help free up resources.
// The iteration's content will be flushed automatically.
// An iteration once closed cannot (yet) be reopened.
i.close();
cout << "Full E/x starts with:\n\t{";
for (size_t col = 0; col < extent[1] && col < 5; ++col)
cout << all_data.get()[col] << ", ";
cout << "...}\n";
/* The files in 'series' are still open until the object is destroyed, on
* which it cleanly flushes and closes all open file handles.
* When running out of scope on return, the 'Series' destructor is called.
* Alternatively, one can call `series.close()` to the same effect as
* calling the destructor, including the release of file handles.
*/
series.close();
return 0;
}
An extended example can be found in examples/6_dump_filebased_series.cpp
.
Python
import openpmd_api as io
if __name__ == "__main__":
series = io.Series("../samples/git-sample/data%T.h5",
io.Access.read_only)
print("Read a Series with openPMD standard version %s" %
series.openPMD)
print("The Series contains {0} iterations:".format(len(series.iterations)))
for i in series.iterations:
print("\t {0}".format(i))
print("")
i = series.iterations[100]
print("Iteration 100 contains {0} meshes:".format(len(i.meshes)))
for m in i.meshes:
print("\t {0}".format(m))
print("")
print("Iteration 100 contains {0} particle species:".format(
len(i.particles)))
for ps in i.particles:
print("\t {0}".format(ps))
print("With records:")
for r in i.particles[ps]:
print("\t {0}".format(r))
# printing a scalar value
electrons = i.particles["electrons"]
charge = electrons["charge"]
series.flush()
print("And the first electron particle has a charge {}"
.format(charge[0]))
print("")
E_x = i.meshes["E"]["x"]
shape = E_x.shape
print("Field E.x has shape {0} and datatype {1}".format(
shape, E_x.dtype))
chunk_data = E_x[1:3, 1:3, 1:2]
# print("Queued the loading of a single chunk from disk, "
# "ready to execute")
series.flush()
print("Chunk has been read from disk\n"
"Read chunk contains:")
print(chunk_data)
# for row in range(2):
# for col in range(2):
# print("\t({0}|{1}|{2})\t{3}".format(
# row + 1, col + 1, 1, chunk_data[row*chunk_extent[1]+col])
# )
# print("")
all_data = E_x.load_chunk()
# The iteration can be closed in order to help free up resources.
# The iteration's content will be flushed automatically.
# An iteration once closed cannot (yet) be reopened.
i.close()
print("Full E/x is of shape {0} and starts with:".format(all_data.shape))
print(all_data[0, 0, :5])
# The files in 'series' are still open until the series is closed, at which
# time it cleanly flushes and closes all open file handles.
# One can close the object explicitly to trigger this.
# Alternatively, this will automatically happen once the garbage collector
# claims (every copy of) the series object.
series.close()
Writing
C++
#include <openPMD/openPMD.hpp>
#include <cstdlib>
#include <iostream>
#include <memory>
#include <numeric>
using std::cout;
using namespace openPMD;
int main(int argc, char *argv[])
{
// user input: size of matrix to write, default 3x3
size_t size = (argc == 2 ? atoi(argv[1]) : 3);
// matrix dataset to write with values 0...size*size-1
std::vector<double> global_data(size * size);
std::iota(global_data.begin(), global_data.end(), 0.);
cout << "Set up a 2D square array (" << size << 'x' << size
<< ") that will be written\n";
// open file for writing
Series series = Series("../samples/3_write_serial.h5", Access::CREATE);
cout << "Created an empty " << series.iterationEncoding() << " Series\n";
// `Series::writeIterations()` and `Series::readIterations()` are
// intentionally restricted APIs that ensure a workflow which also works
// in streaming setups, e.g. an iteration cannot be opened again once
// it has been closed.
// `Series::iterations` can be directly accessed in random-access workflows.
Mesh rho = series.writeIterations()[1].meshes["rho"];
cout << "Created a scalar mesh Record with all required openPMD "
"attributes\n";
Datatype datatype = determineDatatype(global_data.data());
Extent extent = {size, size};
Dataset dataset = Dataset(datatype, extent);
cout << "Created a Dataset of size " << dataset.extent[0] << 'x'
<< dataset.extent[1] << " and Datatype " << dataset.dtype << '\n';
rho.resetDataset(dataset);
cout << "Set the dataset properties for the scalar field rho in iteration "
"1\n";
series.flush();
cout << "File structure and required attributes have been written\n";
Offset offset = {0, 0};
rho.storeChunk(global_data, offset, extent);
cout << "Stored the whole Dataset contents as a single chunk, "
"ready to write content\n";
// The iteration can be closed in order to help free up resources.
// The iteration's content will be flushed automatically.
// An iteration once closed cannot (yet) be reopened.
series.writeIterations()[1].close();
cout << "Dataset content has been fully written\n";
/* The files in 'series' are still open until the object is destroyed, on
* which it cleanly flushes and closes all open file handles.
* When running out of scope on return, the 'Series' destructor is called.
* Alternatively, one can call `series.close()` to the same effect as
* calling the destructor, including the release of file handles.
*/
series.close();
return 0;
}
An extended example can be found in examples/7_extended_write_serial.cpp
.
Python
import numpy as np
import openpmd_api as io
if __name__ == "__main__":
# user input: size of matrix to write, default 3x3
size = 3
# matrix dataset to write with values 0...size*size-1
data = np.arange(size*size, dtype=np.double).reshape(3, 3)
print("Set up a 2D square array ({0}x{1}) that will be written".format(
size, size))
# open file for writing
series = io.Series(
"../samples/3_write_serial_py.h5",
io.Access.create
)
print("Created an empty {0} Series".format(series.iteration_encoding))
print(len(series.iterations))
# `Series.write_iterations()` and `Series.read_iterations()` are
# intentionally restricted APIs that ensure a workflow which also works
# in streaming setups, e.g. an iteration cannot be opened again once
# it has been closed.
# `Series.iterations` can be directly accessed in random-access workflows.
rho = series.write_iterations()[1]. \
meshes["rho"]
dataset = io.Dataset(data.dtype, data.shape)
print("Created a Dataset of size {0}x{1} and Datatype {2}".format(
dataset.extent[0], dataset.extent[1], dataset.dtype))
rho.reset_dataset(dataset)
print("Set the dataset properties for the scalar field rho in iteration 1")
series.flush()
print("File structure has been written")
rho[()] = data
print("Stored the whole Dataset contents as a single chunk, " +
"ready to write content")
# The iteration can be closed in order to help free up resources.
# The iteration's content will be flushed automatically.
# An iteration once closed cannot (yet) be reopened.
series.write_iterations()[1].close()
print("Dataset content has been fully written")
# The files in 'series' are still open until the series is closed, at which
# time it cleanly flushes and closes all open file handles.
# One can close the object explicitly to trigger this.
# Alternatively, this will automatically happen once the garbage collector
# claims (every copy of) the series object.
series.close()