added: ability to put terminal log in HDF5

This commit is contained in:
Arne Morten Kvarving 2020-04-28 10:12:47 +02:00 committed by Kjetil Andre Johannessen
parent bdf5e31c3e
commit 0b63ec3e93
10 changed files with 123 additions and 44 deletions

View File

@ -35,6 +35,7 @@ char** IFEM::argv = nullptr;
SIMoptions IFEM::cmdOptions;
ControlFIFO* IFEM::fifo = nullptr;
utl::LogStream IFEM::cout(std::cout);
std::shared_ptr<std::ostringstream> IFEM::memoryLog;
int IFEM::Init (int arg_c, char** arg_v, const char* title)
@ -52,6 +53,8 @@ int IFEM::Init (int arg_c, char** arg_v, const char* title)
cmdOptions.parseOldOptions(argc,argv,i);
cout.setPIDs(0,linalg.myPid);
memoryLog = std::make_shared<std::ostringstream>();
cout.addExtraLog(memoryLog);
if (linalg.myPid != 0 || argc < 2)
return linalg.myPid;
@ -59,9 +62,9 @@ int IFEM::Init (int arg_c, char** arg_v, const char* title)
if (title)
{
int i, nchar = 13 + strlen(title);
std::cout <<"\n >>> IFEM "<< title <<" <<<\n ";
IFEM::cout <<"\n >>> IFEM "<< title <<" <<<\n ";
for (i = 0; i < nchar; i++) std::cout <<'=';
std::cout <<"\n\n Executing command:\n";
IFEM::cout <<"\n\n Executing command:\n";
#ifdef HAVE_MPI
int nProc = 1;
#ifdef HAS_PETSC
@ -70,91 +73,91 @@ int IFEM::Init (int arg_c, char** arg_v, const char* title)
MPI_Comm_size(MPI_COMM_WORLD,&nProc);
#endif
if (nProc > 1)
std::cout <<" mpirun -np "<< nProc;
IFEM::cout <<" mpirun -np "<< nProc;
#endif
for (i = 0; i < argc; i++) IFEM::cout <<" "<< argv[i];
std::cout << std::endl;
IFEM::cout << std::endl;
}
std::cout <<"\n ===== IFEM v"<< IFEM_VERSION_MAJOR <<"."
<< IFEM_VERSION_MINOR <<"."
<< IFEM_VERSION_PATCH <<" initialized =====";
IFEM::cout <<"\n ===== IFEM v"<< IFEM_VERSION_MAJOR <<"."
<< IFEM_VERSION_MINOR <<"."
<< IFEM_VERSION_PATCH <<" initialized =====";
std::cout <<"\n HDF5 support: ";
IFEM::cout <<"\n HDF5 support: ";
#if HAS_HDF5
std::cout <<"enabled";
IFEM::cout <<"enabled";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n LR spline support: ";
IFEM::cout <<"\n LR spline support: ";
#if HAS_LRSPLINE
std::cout <<"enabled";
IFEM::cout <<"enabled";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n OpenMP support: ";
IFEM::cout <<"\n OpenMP support: ";
#if USE_OPENMP
std::cout <<"enabled";
IFEM::cout <<"enabled";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n MPI support: ";
IFEM::cout <<"\n MPI support: ";
#ifdef HAVE_MPI
std::cout <<"enabled";
IFEM::cout <<"enabled";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n PETSc support: ";
IFEM::cout <<"\n PETSc support: ";
#if HAS_PETSC
std::cout <<"enabled (v"<< PETSC_VERSION_MAJOR <<"."
<< PETSC_VERSION_MINOR <<"."
<< PETSC_VERSION_SUBMINOR <<")";
IFEM::cout <<"enabled (v"<< PETSC_VERSION_MAJOR <<"."
<< PETSC_VERSION_MINOR <<"."
<< PETSC_VERSION_SUBMINOR <<")";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n SuperLU support: ";
IFEM::cout <<"\n SuperLU support: ";
#if HAS_SUPERLU
std::cout <<"enabled (serial)";
IFEM::cout <<"enabled (serial)";
#elif HAS_SUPERLU_MT
std::cout <<"enabled (multi-threaded)";
IFEM::cout <<"enabled (multi-threaded)";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n UMFPack support: ";
IFEM::cout <<"\n UMFPack support: ";
#if HAS_UMFPACK
std::cout <<"enabled (v" << UMFPACK_MAIN_VERSION <<"."
<< UMFPACK_SUB_VERSION <<"."
<< UMFPACK_SUBSUB_VERSION << ")";
IFEM::cout <<"enabled (v" << UMFPACK_MAIN_VERSION <<"."
<< UMFPACK_SUB_VERSION <<"."
<< UMFPACK_SUBSUB_VERSION << ")";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n ISTL support: ";
IFEM::cout <<"\n ISTL support: ";
#if HAS_ISTL
std::cout <<"enabled (v"<< ISTL_VERSION <<")";
IFEM::cout <<"enabled (v"<< ISTL_VERSION <<")";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
std::cout <<"\n VTF support: ";
IFEM::cout <<"\n VTF support: ";
#if HAS_VTFAPI == 2
std::cout <<"enabled (v2)";
IFEM::cout <<"enabled (v2)";
#elif HAS_VTFAPI == 1
std::cout <<"enabled (v1)";
IFEM::cout <<"enabled (v1)";
#else
std::cout <<"disabled";
IFEM::cout <<"disabled";
#endif
if (enableController)
{
fifo = new ControlFIFO();
if (fifo->open())
std::cout <<"\n External controller enabled";
IFEM::cout <<"\n External controller enabled";
else
{
delete fifo;
@ -162,7 +165,7 @@ int IFEM::Init (int arg_c, char** arg_v, const char* title)
}
}
std::cout << std::endl;
IFEM::cout << std::endl;
return 0;
}

View File

@ -61,6 +61,8 @@ public:
static utl::LogStream cout; //!< Combined standard out for parallel processes
static std::shared_ptr<std::ostringstream> memoryLog; //!< Logging data stored in memory
private:
static SIMoptions cmdOptions; //!< General input options
static int argc; //!< The number of command-line arguments

View File

@ -44,7 +44,7 @@ SIMoptions::SIMoptions ()
format = -1;
saveInc = 1;
dtSave = 0.0;
pSolOnly = saveNorms = false;
pSolOnly = saveNorms = saveLog = false;
restartInc = 0;
restartStep = -1;
@ -182,6 +182,9 @@ bool SIMoptions::parseOutputTag (const TiXmlElement* elem)
else if (!strcasecmp(elem->Value(),"saveNorms"))
saveNorms = true;
else if (!strcasecmp(elem->Value(),"saveLog"))
saveLog = true;
else if (!strcasecmp(elem->Value(),"projection")) {
std::string type;
if (utl::getAttribute(elem,"type",type))
@ -248,7 +251,13 @@ bool SIMoptions::parseConsoleTag (const TiXmlElement* elem)
bool SIMoptions::dumpHDF5 (const char* defaultName)
{
if (hdf5.empty()) return false;
if (!saveLog) {
IFEM::cout.removeExtraLog(IFEM::memoryLog);
IFEM::memoryLog.reset();
}
if (hdf5.empty())
return false;
if (hdf5 == "(default)") {
hdf5 = defaultName;

View File

@ -90,6 +90,7 @@ public:
double dtSave; //!< Time interval between each result output
bool pSolOnly; //!< If \e true, don't save secondary solution variables
bool saveNorms;//!< If \e true, save element norms
bool saveLog; //!< If \e true, export the log
std::string hdf5; //!< Prefix for HDF5-file
std::string vtf; //!< Prefix for VTF-file

View File

@ -12,6 +12,7 @@
//==============================================================================
#include "DataExporter.h"
#include "IFEM.h"
#include "Utilities.h"
#include "Profiler.h"
#include "ProcessAdm.h"
@ -35,6 +36,13 @@ DataWriter::DataWriter (const std::string& name,
DataExporter::~DataExporter ()
{
if (IFEM::memoryLog)
for (DataWriter* writer : m_writers) {
writer->openFile(0);
writer->writeLog(IFEM::memoryLog->str(), "cout");
writer->closeFile(0);
}
if (m_delete)
for (DataWriter* writer : m_writers)
delete writer;

View File

@ -16,12 +16,16 @@
#include "ControlFIFO.h"
#include <map>
#include <memory>
#include <string>
#include <vector>
class DataWriter;
class ProcessAdm;
class TimeStep;
namespace utl {
class LogStream;
}
/*!
@ -227,6 +231,10 @@ public:
virtual bool writeTimeInfo(int level, int interval,
const TimeStep& tp) = 0;
//! \brief Write a log to output file.
//! \param name Name of log
virtual bool writeLog(const std::string& data, const std::string& name) = 0;
//! \brief Sets the prefices used for norm output.
void setNormPrefixes(const std::vector<std::string>& prefix) { m_prefix = prefix; }

View File

@ -727,3 +727,34 @@ void HDF5Writer::writeNodalForces(int level, const DataEntry& entry)
std::cout << "HDF5Writer: compiled without HDF5 support, no data written" << std::endl;
#endif
}
bool HDF5Writer::writeLog (const std::string& data, const std::string& name)
{
#ifdef HAS_HDF5
hid_t group;
if (checkGroupExistence(m_file,"/log"))
group = H5Gopen2(m_file,"/log",H5P_DEFAULT);
else
group = H5Gcreate2(m_file,"/log",0,H5P_DEFAULT,H5P_DEFAULT);
int rank = 0;
int size = 1;
#ifdef HAVE_MPI
MPI_Comm_rank(*m_adm.getCommunicator(), &rank);
MPI_Comm_size(*m_adm.getCommunicator(), &size);
#endif
for (int i = 0; i < size; ++i) {
if (rank == i)
writeArray(group, name.c_str(), i, data.size(), data.data(), H5T_NATIVE_CHAR);
else {
char dummy;
writeArray(group, name.c_str(), i, 0, &dummy, H5T_NATIVE_CHAR);
}
}
H5Gclose(group);
return true;
#endif
}

View File

@ -92,6 +92,11 @@ public:
//! \param[in] tp The current time stepping info
virtual bool writeTimeInfo(int level, int interval, const TimeStep& tp);
//! \brief Write a log to output file.
//! \param data Text to write
//! \param name Name of log
virtual bool writeLog(const std::string& data, const std::string& name);
#ifdef HAS_HDF5
protected:
//! \brief Internal helper function writing a data array to file.

View File

@ -13,6 +13,8 @@
#include "LogStream.h"
#include <algorithm>
utl::LogStream::LogStream(std::ostream& out, int ppid, int mypid) :
m_out(&out), m_ppid(ppid), m_pid(mypid)
@ -59,6 +61,14 @@ void utl::LogStream::addExtraLog(std::shared_ptr<std::ostream> extra, bool clear
}
void utl::LogStream::removeExtraLog(std::shared_ptr<std::ostream> extra)
{
auto it = std::find(m_extra.begin(), m_extra.end(), extra);
if (it != m_extra.end())
m_extra.erase(it);
}
int utl::LogStream::precision(int streamsize)
{
int result = streamsize;

View File

@ -45,6 +45,8 @@ public:
void addExtraLog(std::ostream* extra, bool clear = false);
//! \brief Adds an extra logging stream.
void addExtraLog(std::shared_ptr<std::ostream> extra, bool clear = false);
//! \brief Drop an extra logging stream.
void removeExtraLog(std::shared_ptr<std::ostream> extra);
//! \brief Write data to stream
template<typename T>