Skip to content

Commit

Permalink
adding inline engine heat transfer example/test
Browse files Browse the repository at this point in the history
  • Loading branch information
caitlinross committed May 17, 2020
1 parent 9f52c39 commit 0cebe96
Show file tree
Hide file tree
Showing 8 changed files with 414 additions and 0 deletions.
1 change: 1 addition & 0 deletions examples/heatTransfer/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -6,3 +6,4 @@
add_subdirectory(write)
add_subdirectory(read)
add_subdirectory(read_fileonly)
add_subdirectory(inline)
38 changes: 38 additions & 0 deletions examples/heatTransfer/heat_inline.xml
Original file line number Diff line number Diff line change
@@ -0,0 +1,38 @@
<?xml version="1.0"?>
<!-- Config XML file fo the
- heatTransfer_inline_adios2
executables in build/bin
When using this XML config, one must launch the reader and writer
application together with a single mpirun command. They will share
the MPI_COMM_WORLD and create two partitions of the processes.
-->

<adios-config>

<!--====================================
Configuration for the Inline Writer and Reader
====================================-->

<io name="inlineIO">
<engine type="Inline">
<parameter key="verbose" value="5"/>
<parameter key="writerID" value="inlinewriter"/>
<parameter key="readerID" value="inlinereader"/>

<!-- Microseconds (default), Milliseconds, Seconds,
Minutes, Hours -->
<parameter key="ProfileUnits" value="Microseconds"/>
</engine>
</io>

<!--=======================================
Configuration for the Reader output
=======================================-->

<io name="readerOutput">
<engine type="File">
</engine>
</io>
</adios-config>
14 changes: 14 additions & 0 deletions examples/heatTransfer/inline/CMakeLists.txt
Original file line number Diff line number Diff line change
@@ -0,0 +1,14 @@
#------------------------------------------------------------------------------#
# Distributed under the OSI-approved Apache License, Version 2.0. See
# accompanying file Copyright.txt for details.
#------------------------------------------------------------------------------#

add_executable(heatTransfer_inline_adios2
main.cpp
InlineIO.cpp
../write/HeatTransfer.cpp
../write/Settings.cpp
)
target_link_libraries(heatTransfer_inline_adios2
adios2::cxx11_mpi MPI::MPI_C ${CMAKE_THREAD_LIBS_INIT}
)
94 changes: 94 additions & 0 deletions examples/heatTransfer/inline/InlineIO.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,94 @@
/*
* Distributed under the OSI-approved Apache License, Version 2.0. See
* accompanying file Copyright.txt for details.
*
* InlineIO.cpp
*
* Created on: May 2020
* Author: Caitlin Ross
*/

#include "InlineIO.h"

InlineIO::InlineIO(const Settings &s, MPI_Comm comm)
{
std::string writerName = "writer";
std::string readerName = "reader";

ad = adios2::ADIOS(s.configfile, comm);
inlineIO = ad.DeclareIO("inlineIO");

auto params = inlineIO.Parameters();
if (params.find("writerID") != params.end())
{
writerName = params["writerID"];
}
if (params.find("readerID") != params.end())
{
readerName = params["readerID"];
}

if (!inlineIO.InConfigFile())
{
inlineIO.SetEngine("Inline");
inlineIO.SetParameters({{"writerID", writerName}, {"readerID", readerName}});
}

// define T as 2D global array
varT = inlineIO.DefineVariable<double>(
"T",
// Global dimensions
{s.gndx, s.gndy},
// starting offset of the local array in the global space
{s.offsx, s.offsy},
// local size, could be defined later using SetSelection()
{s.ndx, s.ndy});

inlineWriter = inlineIO.Open(writerName, adios2::Mode::Write, comm);
inlineReader = inlineIO.Open(readerName, adios2::Mode::Read, comm);

// Promise that we are not going to change the variable sizes nor add new
// variables
inlineWriter.LockWriterDefinitions();
inlineReader.LockReaderSelections();
}

InlineIO::~InlineIO()
{
inlineWriter.Close();
inlineReader.Close();
}

void InlineIO::write(const HeatTransfer& ht)
{
inlineWriter.BeginStep();
v = ht.data_noghost();
inlineWriter.Put<double>(varT, v.data());
inlineWriter.EndStep();
}

const double* InlineIO::read(bool firstStep)
{
adios2::StepStatus status =
inlineReader.BeginStep(adios2::StepMode::Read);
if (status != adios2::StepStatus::OK)
{
return nullptr;
}

auto blocksInfo = inlineReader.BlocksInfo(varT, inlineReader.CurrentStep());
// in this example we're only expecting one block
if (blocksInfo.size() != 1)
{
throw std::runtime_error("InlineIO::read found incorrect number of blocks");
}

auto& info = blocksInfo[0];
varT.SetBlockSelection(info.BlockID);
inlineReader.Get(varT, info);
inlineReader.EndStep();

// now we can get the pointer with info.Data()
return info.Data();
}

43 changes: 43 additions & 0 deletions examples/heatTransfer/inline/InlineIO.h
Original file line number Diff line number Diff line change
@@ -0,0 +1,43 @@
/*
* Distributed under the OSI-approved Apache License, Version 2.0. See
* accompanying file Copyright.txt for details.
*
* InlineIO.h
*
* Created on: May 2020
* Author: Caitlin Ross
*/

#ifndef InlineIO_H_
#define InlineIO_H_

#include "../write/HeatTransfer.h"
#include "../write/Settings.h"

#include <mpi.h>

#include "adios2.h"

class InlineIO
{
public:
InlineIO(const Settings &s, MPI_Comm comm);
~InlineIO();
void write(const HeatTransfer &ht);
const double* read(bool firstStep);

private:
adios2::ADIOS ad;
adios2::IO inlineIO;
adios2::Engine inlineWriter;
adios2::Engine inlineReader;

adios2::Variable<double> varT;

// need this so data on the write side doesn't go out
// of scope before the reader can use it
std::vector<double> v;

};

#endif /* InlineIO_H_ */
191 changes: 191 additions & 0 deletions examples/heatTransfer/inline/main.cpp
Original file line number Diff line number Diff line change
@@ -0,0 +1,191 @@
/*
* Distributed under the OSI-approved Apache License, Version 2.0. See
* accompanying file Copyright.txt for details.
*
* main.cpp
*
* Recreates heat_transfer.f90 (Fortran) ADIOS tutorial example in C++
* This version is adapted from heatTransfer/write/main.cpp and
* heatTransfer/read/heatRead.cpp for use with the inline engine,
* which requires the writer and reader to be created in the same process
* and the same ADIOS IO object.
*
* Created on: May 2020
* Author: Norbert Podhorszki
* Caitlin Ross
*
*/
#include <mpi.h>

#include "adios2.h"

#include <iostream>
#include <memory>
#include <stdexcept>
#include <string>

#include "InlineIO.h"
#include "../write/HeatTransfer.h"
#include "../write/Settings.h"

void printUsage()
{
std::cout << "Usage: heatTransfer config output N M nx ny steps "
"iterations\n"
<< " config: XML config file to use\n"
<< " output: name of output data file/stream\n"
<< " N: number of processes in X dimension\n"
<< " M: number of processes in Y dimension\n"
<< " nx: local array size in X dimension per processor\n"
<< " ny: local array size in Y dimension per processor\n"
<< " steps: the total number of steps to output\n"
<< " iterations: one step consist of this many iterations\n\n";
}

void Compute(const double* Tin, std::vector<double> &Tout,
std::vector<double> &dT, bool firstStep)
{
/* Compute dT and
* copy Tin into Tout as it will be used for calculating dT in the
* next step
*/
if (firstStep)
{
for (size_t i = 0; i < dT.size(); ++i)
{
dT[i] = 0;
Tout[i] = Tin[i];
}
}
else
{
for (size_t i = 0; i < dT.size(); ++i)
{
dT[i] = Tout[i] - Tin[i];
Tout[i] = Tin[i];
}
}
}

std::vector<double> Tout;
std::vector<double> dT;
adios2::Variable<double> vTout;
adios2::Variable<double> vdT;
adios2::IO outIO;
adios2::ADIOS ad;
adios2::Engine writer;

void writeOutput()
{
writer.BeginStep();
if (vTout)
writer.Put<double>(vTout, Tout.data());
if (vdT)
writer.Put<double>(vdT, dT.data());
writer.EndStep();
}

void setupOutputIO(const Settings& s, MPI_Comm comm)
{
ad = adios2::ADIOS(s.configfile, comm);
outIO = ad.DeclareIO("readerOutput");
Tout.resize(s.ndx * s.ndy);
dT.resize(s.ndx * s.ndy);

/* Create output variables and open output stream */
// For inline engine, there's no exchange of data between processes,
// so the shape of variables to be written out for validation
// is the same as the writer
vTout = outIO.DefineVariable<double>(
"T", {s.gndx, s.gndy}, {s.offsx, s.offsy}, {s.ndx, s.ndy});
vdT = outIO.DefineVariable<double>(
"dT", {s.gndx, s.gndy}, {s.offsx, s.offsy}, {s.ndx, s.ndy});
writer = outIO.Open(s.outputfile, adios2::Mode::Write,
comm);
}

int main(int argc, char *argv[])
{
MPI_Init(&argc, &argv);

/* When writer and reader is launched together with a single mpirun command,
the world comm spans all applications. We have to split and create the
local 'world' communicator mpiHeatTransferComm for the writer only.
When writer and reader is launched separately, the mpiHeatTransferComm
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);

const unsigned int color = 1;
MPI_Comm mpiHeatTransferComm;
MPI_Comm_split(MPI_COMM_WORLD, color, wrank, &mpiHeatTransferComm);

int rank, nproc;
MPI_Comm_rank(mpiHeatTransferComm, &rank);
MPI_Comm_size(mpiHeatTransferComm, &nproc);

try
{
double timeStart = MPI_Wtime();
Settings settings(argc, argv, rank, nproc);
HeatTransfer ht(settings);
InlineIO io(settings, mpiHeatTransferComm);
setupOutputIO(settings, mpiHeatTransferComm);

ht.init(false);
// ht.printT("Initialized T:", mpiHeatTransferComm);
ht.heatEdges();
ht.exchange(mpiHeatTransferComm);
// ht.printT("Heated T:", mpiHeatTransferComm);

io.write(ht);
const double* Tin = io.read(true);
Compute(Tin, Tout, dT, true);
writeOutput();

for (unsigned int t = 1; t <= settings.steps; ++t)
{
if (rank == 0)
std::cout << "Step " << t << ":\n";
for (unsigned int iter = 1; iter <= settings.iterations; ++iter)
{
ht.iterate();
ht.exchange(mpiHeatTransferComm);
ht.heatEdges();
}

io.write(ht);
Tin = io.read(false);
Compute(Tin, Tout, dT, false);
writeOutput();
}
MPI_Barrier(mpiHeatTransferComm);

double timeEnd = MPI_Wtime();
if (rank == 0)
std::cout << "Total runtime = " << timeEnd - timeStart << "s\n";

writer.Close();
}
catch (std::invalid_argument &e) // command-line argument errors
{
std::cout << e.what() << std::endl;
printUsage();
}
catch (std::ios_base::failure &e) // I/O failure (e.g. file not found)
{
std::cout << "I/O base exception caught\n";
std::cout << e.what() << std::endl;
}
catch (std::exception &e) // All other exceptions
{
std::cout << "Exception caught\n";
std::cout << e.what() << std::endl;
}

MPI_Finalize();
return 0;
}
2 changes: 2 additions & 0 deletions testing/examples/heatTransfer/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -7,6 +7,8 @@ include(${CMAKE_CURRENT_SOURCE_DIR}/TestBPFileMx1.cmake)
include(${CMAKE_CURRENT_SOURCE_DIR}/TestBPFileMxM.cmake)
include(${CMAKE_CURRENT_SOURCE_DIR}/TestBPFileMxN.cmake)

include(${CMAKE_CURRENT_SOURCE_DIR}/TestInlineMxM.cmake)

if(ADIOS2_HAVE_ZFP)
#include(${CMAKE_CURRENT_SOURCE_DIR}/TestBPFileMxM_zfp.cmake)
#include(${CMAKE_CURRENT_SOURCE_DIR}/TestBPFileMxN_zfp.cmake)
Expand Down
Loading

0 comments on commit 0cebe96

Please sign in to comment.