Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add HDF5 input and output files to eigensolver miniapp #973

Merged
merged 11 commits into from
Sep 13, 2023
55 changes: 49 additions & 6 deletions miniapp/miniapp_eigensolver.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -27,6 +27,7 @@
#include <dlaf/eigensolver/internal/get_band_size.h>
#include <dlaf/init.h>
#include <dlaf/matrix/copy.h>
#include <dlaf/matrix/hdf5.h>
#include <dlaf/matrix/matrix.h>
#include <dlaf/matrix/matrix_mirror.h>
#include <dlaf/miniapp/dispatch.h>
Expand All @@ -50,6 +51,10 @@ using dlaf::common::Ordering;
using dlaf::matrix::MatrixMirror;
using pika::this_thread::experimental::sync_wait;

#ifdef DLAF_WITH_HDF5
using dlaf::matrix::internal::FileHDF5;
#endif

/// Check results of the eigensolver
template <typename T>
void checkEigensolver(CommunicatorGrid comm_grid, blas::Uplo uplo, Matrix<const T, Device::CPU>& A,
Expand All @@ -60,12 +65,31 @@ struct Options
SizeType m;
SizeType mb;
blas::Uplo uplo;
#ifdef DLAF_WITH_HDF5
std::optional<dlaf::matrix::internal::FileHDF5> input_file;
std::optional<dlaf::matrix::internal::FileHDF5> output_file;
#endif

Options(const pika::program_options::variables_map& vm)
: MiniappOptions(vm), m(vm["matrix-size"].as<SizeType>()), mb(vm["block-size"].as<SizeType>()),
uplo(dlaf::miniapp::parseUplo(vm["uplo"].as<std::string>())) {
DLAF_ASSERT(m > 0, m);
DLAF_ASSERT(mb > 0, mb);

#ifdef DLAF_WITH_HDF5
if (vm.count("input-file") == 1) {
input_file = FileHDF5(vm["input-file"].as<std::string>(), FileHDF5::FileMode::readonly);

if (!vm["matrix-size"].defaulted()) {
std::cerr << "Warning! "
"Specified matrix size will be ignored because an input file has been specified."
<< std::endl;
}
}
if (vm.count("output-file") == 1) {
output_file = FileHDF5(vm["output-file"].as<std::string>(), FileHDF5::FileMode::readwrite);
}
#endif
}

Options(Options&&) = default;
Expand All @@ -87,19 +111,25 @@ struct EigensolverMiniapp {
Communicator world(MPI_COMM_WORLD);
CommunicatorGrid comm_grid(world, opts.grid_rows, opts.grid_cols, Ordering::ColumnMajor);

// Allocate memory for the matrix
GlobalElementSize matrix_size(opts.m, opts.m);
TileElementSize block_size(opts.mb, opts.mb);

ConstHostMatrixType matrix_ref = [matrix_size, block_size, comm_grid]() {
ConstHostMatrixType matrix_ref = [comm_grid, &opts]() {
#ifdef DLAF_WITH_HDF5
if (opts.input_file) {
Matrix<T, Device::CPU> input = opts.input_file->read<T>("/a", {opts.mb, opts.mb});
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Suggested change
Matrix<T, Device::CPU> input = opts.input_file->read<T>("/a", {opts.mb, opts.mb});
Matrix<T, Device::CPU> input = opts.input_file->read<T>("/input", {opts.mb, opts.mb});

?

return input;
}
#endif
using dlaf::matrix::util::set_random_hermitian;

HostMatrixType hermitian(matrix_size, block_size, comm_grid);
HostMatrixType hermitian(GlobalElementSize(opts.m, opts.m), TileElementSize(opts.mb, opts.mb),
comm_grid);
set_random_hermitian(hermitian);

return hermitian;
}();

auto matrix_size = matrix_ref.size();
auto block_size = matrix_ref.blockSize();

for (int64_t run_index = -opts.nwarmups; run_index < opts.nruns; ++run_index) {
if (0 == world.rank() && run_index >= 0)
std::cout << "[" << run_index << "]" << std::endl;
Expand Down Expand Up @@ -127,6 +157,15 @@ struct EigensolverMiniapp {
DLAF_MPI_CHECK_ERROR(MPI_Barrier(world));
double elapsed_time = timeit.elapsed();

#ifdef DLAF_WITH_HDF5
if (run_index == opts.nruns - 1) {
if (opts.output_file) {
opts.output_file->write<BaseType<T>>(eigenvalues, "/evals");
opts.output_file->write<T>(eigenvectors, "/evecs");
}
}
#endif

matrix.reset();

// print benchmark results
Expand Down Expand Up @@ -176,6 +215,10 @@ int main(int argc, char** argv) {
desc_commandline.add_options()
("matrix-size", value<SizeType>() ->default_value(4096), "Matrix size")
("block-size", value<SizeType>() ->default_value( 256), "Block cyclic distribution size")
#ifdef DLAF_WITH_HDF5
("input-file", value<std::string>() , "Load matrix from given HDF5 file")
("output-file", value<std::string>() , "Save eigenvector and eigenvalue to given HDF5 file")
RMeli marked this conversation as resolved.
Show resolved Hide resolved
#endif
;
// clang-format on
dlaf::miniapp::addUploOption(desc_commandline);
Expand Down
Loading