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 history file output #366

Merged
merged 4 commits into from
Sep 8, 2023
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
10 changes: 10 additions & 0 deletions src/AdvectionSimulation.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -82,6 +82,9 @@ template <typename problem_t> class AdvectionSimulation : public AMRSimulation<p
// compute derived variables
void ComputeDerivedVar(int lev, std::string const &dname, amrex::MultiFab &mf, int ncomp) const override;

// compute statistics
auto ComputeStatistics() -> std::map<std::string, amrex::Real> override;

void FixupState(int lev) override;

// tag cells for refinement
Expand Down Expand Up @@ -152,6 +155,13 @@ template <typename problem_t> void AdvectionSimulation<problem_t>::ComputeDerive
// user should implement
}

template <typename problem_t>
auto AdvectionSimulation<problem_t>::ComputeStatistics() -> std::map<std::string, amrex::Real>
{
// user should implement
return std::map<std::string, amrex::Real>{};
}

template <typename problem_t> void AdvectionSimulation<problem_t>::ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real /*time*/, int /*ngrow*/)
{
// tag cells for refinement -- implement in problem generator
Expand Down
11 changes: 11 additions & 0 deletions src/RadhydroSimulation.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -191,6 +191,9 @@ template <typename problem_t> class RadhydroSimulation : public AMRSimulation<pr
// compute derived variables
void ComputeDerivedVar(int lev, std::string const &dname, amrex::MultiFab &mf, int ncomp) const override;

// compute statistics
auto ComputeStatistics() -> std::map<std::string, amrex::Real> override;

// fix-up states
void FixupState(int level) override;

Expand Down Expand Up @@ -496,6 +499,14 @@ void RadhydroSimulation<problem_t>::ComputeDerivedVar(int lev, std::string const
// compute derived variables and save in 'mf' -- user should implement
}

template <typename problem_t>
auto RadhydroSimulation<problem_t>::ComputeStatistics() -> std::map<std::string, amrex::Real>
{
// compute statistics and return a std::map<std::string, amrex::Real> -- user should implement
// IMPORTANT: the user is responsible for performing any necessary MPI reductions before returning
markkrumholz marked this conversation as resolved.
Show resolved Hide resolved
return std::map<std::string, amrex::Real>{};
}

template <typename problem_t> void RadhydroSimulation<problem_t>::ErrorEst(int lev, amrex::TagBoxArray &tags, amrex::Real /*time*/, int /*ngrow*/)
{
// tag cells for refinement -- user should implement
Expand Down
69 changes: 69 additions & 0 deletions src/simulation.hpp
Original file line number Diff line number Diff line change
Expand Up @@ -114,6 +114,7 @@ template <typename problem_t> class AMRSimulation : public amrex::AmrCore
amrex::Long maxWalltime_ = 0; // default: no limit
int ascentInterval_ = -1; // -1 == no in-situ renders with Ascent
int plotfileInterval_ = -1; // -1 == no output
int statisticsInterval_ = -1; // -1 == no output
amrex::Real plotTimeInterval_ = -1.0; // time interval for plt file
amrex::Real checkpointTimeInterval_ = -1.0; // time interval for checkpoints
int checkpointInterval_ = -1; // -1 == no output
Expand Down Expand Up @@ -158,6 +159,9 @@ template <typename problem_t> class AMRSimulation : public amrex::AmrCore
// compute derived variables
virtual void ComputeDerivedVar(int lev, std::string const &dname, amrex::MultiFab &mf, int ncomp) const = 0;

// compute statistics
virtual auto ComputeStatistics() -> std::map<std::string, amrex::Real> = 0;

// fix-up any unphysical states created by AMR operations
// (e.g., caused by the flux register or from interpolation)
virtual void FixupState(int level) = 0;
Expand Down Expand Up @@ -222,6 +226,7 @@ template <typename problem_t> class AMRSimulation : public amrex::AmrCore
[[nodiscard]] auto PlotFileMFAtLevel(int lev) -> amrex::MultiFab;
void WriteMetadataFile(std::string const &plotfilename) const;
void ReadMetadataFile(std::string const &chkfilename);
void WriteStatisticsFile();
void WritePlotFile();
void WriteCheckpointFile() const;
void SetLastCheckpointSymlink(std::string const &checkpointname) const;
Expand Down Expand Up @@ -266,6 +271,7 @@ template <typename problem_t> class AMRSimulation : public amrex::AmrCore
/// output parameters
std::string plot_file{"plt"}; // plotfile prefix
std::string chk_file{"chk"}; // checkpoint prefix
std::string stats_file{"history.txt"}; // statistics filename
/// input parameters (if >= 0 we restart from a checkpoint)
std::string restart_chkfile;

Expand Down Expand Up @@ -433,6 +439,9 @@ template <typename problem_t> void AMRSimulation<problem_t>::readParameters()
// Default output interval
pp.query("plotfile_interval", plotfileInterval_);

// Default statistics interval
pp.query("statistics_interval", statisticsInterval_);

// Default Time interval
pp.query("plottime_interval", plotTimeInterval_);

Expand Down Expand Up @@ -524,6 +533,10 @@ template <typename problem_t> void AMRSimulation<problem_t>::setInitialCondition
WritePlotFile();
}

if (statisticsInterval_ > 0) {
WriteStatisticsFile();
}

// ensure that there are enough boxes per MPI rank
PerformanceHints();
}
Expand Down Expand Up @@ -662,6 +675,7 @@ template <typename problem_t> void AMRSimulation<problem_t>::evolve()
#ifdef AMREX_USE_ASCENT
int last_ascent_step = 0;
#endif
int last_statistics_step = 0;
int last_plot_file_step = 0;
double next_plot_file_time = plotTimeInterval_;
double next_chk_file_time = checkpointTimeInterval_;
Expand Down Expand Up @@ -715,6 +729,11 @@ template <typename problem_t> void AMRSimulation<problem_t>::evolve()
}
#endif

if (statisticsInterval_ > 0 && (step + 1) % statisticsInterval_ == 0) {
last_statistics_step = step + 1;
WriteStatisticsFile();
}

if (plotfileInterval_ > 0 && (step + 1) % plotfileInterval_ == 0) {
last_plot_file_step = step + 1;
WritePlotFile();
Expand Down Expand Up @@ -786,6 +805,11 @@ template <typename problem_t> void AMRSimulation<problem_t>::evolve()
WritePlotFile();
}

// write final statistics
if (statisticsInterval_ > 0 && istep[0] > last_statistics_step) {
WriteStatisticsFile();
}

#ifdef AMREX_USE_ASCENT
// close Ascent
ascent_.close();
Expand Down Expand Up @@ -1754,6 +1778,51 @@ template <typename problem_t> void AMRSimulation<problem_t>::WritePlotFile()
#endif
}

template <typename problem_t>
void AMRSimulation<problem_t>::WriteStatisticsFile() {
// append to statistics file
static bool isHeaderWritten = false;

// compute statistics
// IMPORTANT: the user is responsible for performing any necessary MPI reductions inside ComputeStatistics
markkrumholz marked this conversation as resolved.
Show resolved Hide resolved
std::map<std::string, amrex::Real> statistics = ComputeStatistics();

// write to file
if (amrex::ParallelDescriptor::IOProcessor()) {
amrex::VisMF::IO_Buffer io_buffer(amrex::VisMF::IO_Buffer_Size);
std::ofstream StatisticsFile;
StatisticsFile.rdbuf()->pubsetbuf(io_buffer.dataPtr(), io_buffer.size());
StatisticsFile.open(stats_file.c_str(), std::ofstream::out |
std::ofstream::app);
if (!StatisticsFile.good()) {
amrex::FileOpenFailed(stats_file);
}

// write header
if (!isHeaderWritten) {
const std::time_t t = std::chrono::system_clock::to_time_t(std::chrono::system_clock::now());
const std::tm now = *std::localtime(&t); // NOLINT(concurrency-mt-unsafe)
StatisticsFile << "## Simulation restarted at: " << std::put_time(&now, "%c %Z") << "\n";
StatisticsFile << "# cycle time ";
for (auto const &[key, value] : statistics) {
StatisticsFile << key << " ";
}
StatisticsFile << "\n";
isHeaderWritten = true;
}

// save statistics to file
StatisticsFile << istep[0] << " "; // cycle
StatisticsFile << tNew_[0] << " "; // time
for (auto const &[key, value] : statistics) {
StatisticsFile << value << " ";
}
StatisticsFile << "\n";

// file closed automatically by destructor
}
}

template <typename problem_t> void AMRSimulation<problem_t>::SetLastCheckpointSymlink(std::string const &checkpointname) const
{
// creates a symlink pointing to the most recent checkpoint
Expand Down