From a7af19d6eac67eda8a5d30f9e64a8300c3c5a760 Mon Sep 17 00:00:00 2001 From: Norbert Podhorszki <pnorbert@ornl.gov> Date: Fri, 1 Dec 2017 09:51:19 -0500 Subject: [PATCH] refactor heatTransfer example: read_fileonly/ contains the file reader examples, read/ contains the streaming version, writer and reader is using heat.xml config file from the root directory --- examples/heatTransfer/CMakeLists.txt | 3 +- examples/heatTransfer/ReadMe.md | 20 ++- examples/heatTransfer/heat.xml | 84 +++++++++ examples/heatTransfer/read/CMakeLists.txt | 16 +- examples/heatTransfer/read/PrintDataStep.h | 63 +++++++ examples/heatTransfer/read/heatRead.cpp | 159 ++++++++++++++++++ .../heatTransfer/read_fileonly/CMakeLists.txt | 19 +++ .../{read => read_fileonly}/PrintData.h | 0 .../{read => read_fileonly}/heatRead_a2h5.cpp | 0 .../heatRead_adios1.cpp | 0 .../heatRead_adios2.cpp | 3 +- 11 files changed, 346 insertions(+), 21 deletions(-) create mode 100644 examples/heatTransfer/heat.xml create mode 100644 examples/heatTransfer/read/PrintDataStep.h create mode 100644 examples/heatTransfer/read/heatRead.cpp create mode 100644 examples/heatTransfer/read_fileonly/CMakeLists.txt rename examples/heatTransfer/{read => read_fileonly}/PrintData.h (100%) rename examples/heatTransfer/{read => read_fileonly}/heatRead_a2h5.cpp (100%) rename examples/heatTransfer/{read => read_fileonly}/heatRead_adios1.cpp (100%) rename examples/heatTransfer/{read => read_fileonly}/heatRead_adios2.cpp (97%) diff --git a/examples/heatTransfer/CMakeLists.txt b/examples/heatTransfer/CMakeLists.txt index f6c932be4..52f9cab1a 100644 --- a/examples/heatTransfer/CMakeLists.txt +++ b/examples/heatTransfer/CMakeLists.txt @@ -4,4 +4,5 @@ #------------------------------------------------------------------------------# add_subdirectory(write) -add_subdirectory(read) \ No newline at end of file +add_subdirectory(read) +add_subdirectory(read_fileonly) diff --git a/examples/heatTransfer/ReadMe.md b/examples/heatTransfer/ReadMe.md index 639eed4e7..eabafd502 100644 --- a/examples/heatTransfer/ReadMe.md +++ b/examples/heatTransfer/ReadMe.md @@ -1,15 +1,21 @@ examples/heatTransfer Test that solves a 2D Poisson equation for temperature in homogeneous media -using finite differences. This examples shows a straight-forward way to hook an -application to the ADIOS2 library for its IO. +using finite differences. This examples shows a straight-forward way to hook +an application to the ADIOS2 library for its IO. -1. read: illustrates an experimental Read API, uses adios1 underneath +1. read: illustrates the Read API that allows running the reader either as -2. write: illustrates the Write API, resulting binaries under build/bin require arguments + * post-mortem to read all output steps + * in situ to read step by step as the writer outputs them + +2. write: illustrates the Write API as well as has implementations of other IO libraries - * adios1 - * adios2 + * adios 1.x * hdf5 - * phdf5 \ No newline at end of file + * phdf5 + +3. read_fileonly: illustrates reading all output steps at once (a single read + statement) into a single contiguous memory block. This approach only works + for post-mortem processing. diff --git a/examples/heatTransfer/heat.xml b/examples/heatTransfer/heat.xml new file mode 100644 index 000000000..56cd16dbb --- /dev/null +++ b/examples/heatTransfer/heat.xml @@ -0,0 +1,84 @@ +<?xml version="1.0"?> +<!-- Config XML file fo the + - heatTransfer_write_adios2 + - heatTransfer_read + executables in build/bin --> + +<adios-config> + + <!--==================================== + Configuration for the Writer + ====================================--> + + <io name="writer"> + <engine type="BPFileWriter"> + + <!-- for vectorized memory operations, make sure your system + enables threads--> + <parameter key="Threads" value="2"/> + + <!-- Microseconds (default), Milliseconds, Seconds, + Minutes, Hours --> + <parameter key="ProfileUnits" value="Microseconds"/> + + <!-- XXKb, XXMb, or XXXGb supported, default=16Kb + (applications might choose an optimal value) --> + <parameter key="InitialBufferSize" value="16Kb"/> + + <!-- XXKb, XXMb, or XXXGb supported, default=Unlimited (until + fails), maximum at each time step + (applications might choose an optimal value) --> + <!-- <parameter key="MaxBufferSize" value="2Gb"/> --> + + <!-- exponential growth factor > 1, default = 1.05 + 1.05 is good for a few large variables, for many small + variables increase the value to 1.5 to 2 + (optimal value is application dependent)--> + <!-- <parameter key="BufferGrowthFactor" value="1.05"/> --> + + </engine> + + <transport type="File"> + + <!-- POSIX, stdio (C FILE*), fstream (C++) --> + <parameter key="Library" value="POSIX"/> + + <!-- For read/write, Microseconds (default), Milliseconds, Seconds, + Minutes, Hours. open/close always in Microseconds --> + <parameter key="ProfileUnits" value="Microseconds"/> + + </transport> + + </io> + + + <!--==================================== + Configuration for the Reader + ====================================--> + + <io name="reader"> + <engine type="BPFileReader"> + + <!-- for vectorized memory operations, make sure your system + enables threads--> + <parameter key="Threads" value="2"/> + + <!-- Microseconds (default), Milliseconds, Seconds, + Minutes, Hours --> + <parameter key="ProfileUnits" value="Microseconds"/> + + </engine> + + <transport type="File"> + + <!-- POSIX, stdio (C FILE*), fstream (C++) --> + <parameter key="Library" value="POSIX"/> + + <!-- For read/write, Microseconds (default), Milliseconds, Seconds, + Minutes, Hours. open/close always in Microseconds --> + <parameter key="ProfileUnits" value="Microseconds"/> + + </transport> + + </io> +</adios-config> diff --git a/examples/heatTransfer/read/CMakeLists.txt b/examples/heatTransfer/read/CMakeLists.txt index 2a4ff60f7..e3e143b83 100644 --- a/examples/heatTransfer/read/CMakeLists.txt +++ b/examples/heatTransfer/read/CMakeLists.txt @@ -4,16 +4,8 @@ #------------------------------------------------------------------------------# if(ADIOS2_HAVE_MPI) - add_executable(heatTransfer_read_adios2 heatRead_adios2.cpp PrintData.h) - target_link_libraries(heatTransfer_read_adios2 adios2 MPI::MPI_C) - - if(ADIOS2_HAVE_HDF5) - add_executable(heatTransfer_read_a2h5 heatRead_a2h5.cpp PrintData.h) - target_link_libraries(heatTransfer_read_a2h5 adios2 MPI::MPI_C) - endif() - - if(ADIOS2_HAVE_ADIOS1) - add_executable(heatTransfer_read_adios1 heatRead_adios1.cpp PrintData.h) - target_link_libraries(heatTransfer_read_adios1 adios1::adios MPI::MPI_C) - endif() + add_executable(heatTransfer_read heatRead.cpp PrintDataStep.h) + target_link_libraries(heatTransfer_read adios2 MPI::MPI_C) + target_compile_definitions(heatTransfer_read PRIVATE + -DDEFAULT_CONFIG=${CMAKE_CURRENT_SOURCE_DIR}/../heat.xml) endif() diff --git a/examples/heatTransfer/read/PrintDataStep.h b/examples/heatTransfer/read/PrintDataStep.h new file mode 100644 index 000000000..17a6a737d --- /dev/null +++ b/examples/heatTransfer/read/PrintDataStep.h @@ -0,0 +1,63 @@ +/* + * Distributed under the OSI-approved Apache License, Version 2.0. See + * accompanying file Copyright.txt for details. + * + * PrintData.h + * + * Created on: Apr 2017 + * Author: Norbert Podhorszki + */ + +#ifndef PRINTDATASTEP_H_ +#define PRINTDATASTEP_H_ + +#include <cstdint> +#include <fstream> +#include <iomanip> +#include <iostream> +#include <string> + +template <class T> +void printDataStep(double *xy, T *size, T *offset, int rank, int step) +{ + std::ofstream myfile; + std::string filename = "data." + std::to_string(rank); + if (step == 0) + { + myfile.open(filename); + } + else + { + myfile.open(filename, std::ios::app); + } + double *data = xy; + uint64_t nelems = size[0] * size[1]; + myfile << "rank=" << rank << " size=" << size[0] << "x" << size[1] + << " offsets=" << offset[0] << ":" << offset[1] << " step=" << step + << std::endl; + + myfile << " time row columns " << offset[1] << "..." + << offset[1] + size[1] - 1 << std::endl; + myfile << " "; + for (int j = 0; j < size[1]; j++) + { + myfile << std::setw(9) << offset[1] + j; + } + myfile << std::endl; + myfile << "------------------------------------------------------------" + "--\n"; + for (int i = 0; i < size[0]; i++) + { + myfile << std::setw(5) << step << std::setw(5) << offset[0] + i; + for (int j = 0; j < size[1]; j++) + { + myfile << std::setw(9) << std::setprecision(2) + << data[i * size[1] + j]; + } + myfile << std::endl; + } + data += nelems; + myfile.close(); +} + +#endif /* PRINTDATASTEP_H_ */ diff --git a/examples/heatTransfer/read/heatRead.cpp b/examples/heatTransfer/read/heatRead.cpp new file mode 100644 index 000000000..83bf75b66 --- /dev/null +++ b/examples/heatTransfer/read/heatRead.cpp @@ -0,0 +1,159 @@ +/* + * Distributed under the OSI-approved Apache License, Version 2.0. See + * accompanying file Copyright.txt for details. + * + * IO_ADIOS2.cpp + * + * Created on: Nov 2017 + * Author: Norbert Podhorszki + * + */ +#include <mpi.h> + +#include "adios2.h" + +#include <cstdint> +#include <iomanip> +#include <iostream> +#include <math.h> +#include <memory> +#include <stdexcept> +#include <string> +#include <vector> + +#include "PrintDataStep.h" + +#define str_helper(X) #X +#define str(X) str_helper(X) +#ifndef DEFAULT_CONFIG +#define DEFAULT_CONFIG "../heat.xml" +#endif +#define DEFAULT_CONFIG_STR str(DEFAULT_CONFIG) + +int main(int argc, char *argv[]) +{ + MPI_Init(&argc, &argv); + + if (argc < 2) + { + std::cout << "Not enough arguments: need an input file\n"; + return 1; + } + const char *inputfile = argv[1]; + + /* World comm spans all applications started with the same aprun command + on a Cray XK6. So we have to split and create the local + 'world' communicator for the reader only. + In normal start-up, the communicator will just equal the MPI_COMM_WORLD. + */ + + int wrank, wnproc; + MPI_Comm_rank(MPI_COMM_WORLD, &wrank); + MPI_Comm_size(MPI_COMM_WORLD, &wnproc); + MPI_Barrier(MPI_COMM_WORLD); + + const unsigned int color = 2; + MPI_Comm mpiReaderComm; + MPI_Comm_split(MPI_COMM_WORLD, color, wrank, &mpiReaderComm); + + int rank, nproc; + MPI_Comm_rank(mpiReaderComm, &rank); + MPI_Comm_size(mpiReaderComm, &nproc); + + adios2::ADIOS ad(std::string(DEFAULT_CONFIG_STR), mpiReaderComm, + adios2::DebugON); + + // Define method for engine creation + // 1. Get method def from config file or define new one + + adios2::IO &bpReaderIO = ad.DeclareIO("reader"); + if (!bpReaderIO.InConfigFile()) + { + // if not defined by user, we can change the default settings + // BPFileWriter is the default engine + bpReaderIO.SetEngine("ADIOS1Reader"); + bpReaderIO.SetParameters({{"num_threads", "2"}}); + + // ISO-POSIX file output is the default transport (called "File") + // Passing parameters to the transport + bpReaderIO.AddTransport("File", {{"verbose", "4"}}); + } + + adios2::Engine &bpReader = + bpReaderIO.Open(inputfile, adios2::Mode::Read, mpiReaderComm); + + unsigned int gndx; + unsigned int gndy; + double *T; + adios2::Dims readsize; + adios2::Dims offset; + adios2::Variable<double> *vT = nullptr; + bool firstStep = true; + int step = 0; + + while (true) + { + adios2::StepStatus status = + bpReader.BeginStep(adios2::StepMode::NextAvailable); + if (status != adios2::StepStatus::OK) + break; + + if (firstStep) + { + adios2::Variable<unsigned int> *vgndx = + bpReaderIO.InquireVariable<unsigned int>("gndx"); + gndx = vgndx->m_Value; + // bpReader.GetSync<unsigned int>("gndx", gndx); + + adios2::Variable<unsigned int> *vgndy = + bpReaderIO.InquireVariable<unsigned int>("gndy"); + gndy = vgndy->m_Value; + // bpReader.GetSync<unsigned int>("gndy", gndy); + + if (rank == 0) + { + std::cout << "gndx = " << gndx << std::endl; + std::cout << "gndy = " << gndy << std::endl; + } + + // 1D decomposition of the columns, which is inefficient for + // reading! + readsize.push_back(gndx); + readsize.push_back(gndy / nproc); + offset.push_back(0LL); + offset.push_back(rank * readsize[1]); + if (rank == nproc - 1) + { + // last process should read all the rest of columns + readsize[1] = gndy - readsize[1] * (nproc - 1); + } + + std::cout << "rank " << rank << " reads " << readsize[1] + << " columns from offset " << offset[1] << std::endl; + + vT = bpReaderIO.InquireVariable<double>("T"); + T = new double[readsize[0] * readsize[1]]; + + // Create a 2D selection for the subset + vT->SetSelection(adios2::Box<adios2::Dims>(offset, readsize)); + firstStep = false; + } + + if (!rank) + { + std::cout << "Processing step " << step << std::endl; + } + // Arrays are read by scheduling one or more of them + // and performing the reads at once + bpReader.GetDeferred<double>(*vT, T); + bpReader.PerformGets(); + + printDataStep(T, readsize.data(), offset.data(), rank, step); + bpReader.EndStep(); + step++; + } + bpReader.Close(); + delete[] T; + MPI_Finalize(); + return 0; +} diff --git a/examples/heatTransfer/read_fileonly/CMakeLists.txt b/examples/heatTransfer/read_fileonly/CMakeLists.txt new file mode 100644 index 000000000..46c5a1e35 --- /dev/null +++ b/examples/heatTransfer/read_fileonly/CMakeLists.txt @@ -0,0 +1,19 @@ +#------------------------------------------------------------------------------# +# Distributed under the OSI-approved Apache License, Version 2.0. See +# accompanying file Copyright.txt for details. +#------------------------------------------------------------------------------# + +if(ADIOS2_HAVE_MPI) + add_executable(heatTransfer_read_fileonly_adios2 heatRead_adios2.cpp PrintData.h) + target_link_libraries(heatTransfer_read_fileonly_adios2 adios2 MPI::MPI_C) + + if(ADIOS2_HAVE_HDF5) + add_executable(heatTransfer_read_fileonly_a2h5 heatRead_a2h5.cpp PrintData.h) + target_link_libraries(heatTransfer_read_fileonly_a2h5 adios2 MPI::MPI_C) + endif() + + if(ADIOS2_HAVE_ADIOS1) + add_executable(heatTransfer_read_fileonly_adios1 heatRead_adios1.cpp PrintData.h) + target_link_libraries(heatTransfer_read_fileonly_adios1 adios1::adios MPI::MPI_C) + endif() +endif() diff --git a/examples/heatTransfer/read/PrintData.h b/examples/heatTransfer/read_fileonly/PrintData.h similarity index 100% rename from examples/heatTransfer/read/PrintData.h rename to examples/heatTransfer/read_fileonly/PrintData.h diff --git a/examples/heatTransfer/read/heatRead_a2h5.cpp b/examples/heatTransfer/read_fileonly/heatRead_a2h5.cpp similarity index 100% rename from examples/heatTransfer/read/heatRead_a2h5.cpp rename to examples/heatTransfer/read_fileonly/heatRead_a2h5.cpp diff --git a/examples/heatTransfer/read/heatRead_adios1.cpp b/examples/heatTransfer/read_fileonly/heatRead_adios1.cpp similarity index 100% rename from examples/heatTransfer/read/heatRead_adios1.cpp rename to examples/heatTransfer/read_fileonly/heatRead_adios1.cpp diff --git a/examples/heatTransfer/read/heatRead_adios2.cpp b/examples/heatTransfer/read_fileonly/heatRead_adios2.cpp similarity index 97% rename from examples/heatTransfer/read/heatRead_adios2.cpp rename to examples/heatTransfer/read_fileonly/heatRead_adios2.cpp index f4a595572..a0833fa10 100644 --- a/examples/heatTransfer/read/heatRead_adios2.cpp +++ b/examples/heatTransfer/read_fileonly/heatRead_adios2.cpp @@ -44,7 +44,7 @@ int main(int argc, char *argv[]) MPI_Comm_rank(mpiReaderComm, &rank); MPI_Comm_size(mpiReaderComm, &nproc); - adios2::ADIOS ad("adios2.xml", mpiReaderComm, adios2::DebugON); + adios2::ADIOS ad(mpiReaderComm, adios2::DebugON); // Define method for engine creation // 1. Get method def from config file or define new one @@ -56,6 +56,7 @@ int main(int argc, char *argv[]) // BPFileWriter is the default engine bpReaderIO.SetEngine("ADIOS1Reader"); bpReaderIO.SetParameters({{"num_threads", "2"}}); + bpReaderIO.SetParameter("OpenAsFile", "true"); // ISO-POSIX file is the default transport // Passing parameters to the transport -- GitLab