1
0
Fork 0
mirror of https://github.com/cosmo-sims/monofonIC.git synced 2024-09-19 17:03:45 +02:00

some beautification of output

This commit is contained in:
Oliver Hahn 2023-06-23 11:38:34 +02:00
parent 736c57c7dd
commit 9ffe2ce991
4 changed files with 46 additions and 10 deletions

View file

@ -293,7 +293,7 @@ namespace cosmology
{
cosmo_param_.set("pnorm", 1.0 / Dplus_target_ / Dplus_target_);
auto sigma8 = this->compute_sigma8();
music::ilog << "Measured sigma_8 for given PS normalisation is " << sigma8 << std::endl;
music::COLUMN_ilog("Measured sigma_8", sigma8 );
}
cosmo_param_.set("sqrtpnorm", std::sqrt(cosmo_param_["pnorm"]));

View file

@ -21,6 +21,8 @@
#include <cstdarg>
#include <fstream>
#include <iostream>
#include <iomanip> // for std::setw
#include <string> // for std::string
namespace music {
@ -148,4 +150,32 @@ extern log_stream wlog;
extern log_stream ilog;
extern log_stream dlog;
// global functions for convenience
inline void COLUMN_ilog(const std::string& text1, const std::string& text2) {
music::ilog << std::setw(32) << std::left << text1 << " : " << text2 << std::endl;
}
inline void COLUMN_ilog(const std::string& text1, bool flag) {
music::ilog << std::setw(32) << std::left << text1 << " : " << (flag?"enabled":"disabled") << std::endl;
}
inline void COLUMN_ilog(const std::string& text1, double value) {
music::ilog << std::setw(32) << std::left << text1 << " : " << value << std::endl;
}
inline void COLUMN_ilog(const std::string& text1, size_t value) {
music::ilog << std::setw(32) << std::left << text1 << " : " << value << std::endl;
}
inline void CODE_PART_msg(const std::string& text) {
music::ilog << "-------------------------------------------------------------------------------" << std::endl;
music::ilog << "|| " << std::setw(48) << std::setfill('.') << text << std::setw(28) << std::right << "||" << std::setfill(' ') << std::endl;
}
inline void CODE_PART_COMPLETE_msg( double wtime ){
music::ilog << "||" << std::setw(65) << std::right << "took : " << std::setw(8) << wtime << "s ||" << std::endl;
music::ilog << "-------------------------------------------------------------------------------" << std::endl;
music::ilog << std::endl;
}
} // namespace music

View file

@ -23,16 +23,17 @@
void memory_report(void)
{
//... report memory usage
size_t curr_mem_high_mark = 0;
local_mem_high_mark = memory::getCurrentRSS();
#if defined(USE_MPI)
MPI_Allreduce(&local_mem_high_mark, &curr_mem_high_mark, 1, MPI_UNSIGNED_LONG_LONG, MPI_MAX, MPI_COMM_WORLD);
#else
curr_mem_high_mark = local_mem_high_mark;
#endif
if( curr_mem_high_mark > 1.1*global_mem_high_mark ){
music::ilog << std::setw(57) << std::setfill(' ') << std::right << "mem high: " << std::setw(8) << curr_mem_high_mark/(1ull<<20) << " MBytes / task" << std::endl;
if (curr_mem_high_mark > 1.1 * global_mem_high_mark) {
music::ilog << std::setw(32) << std::setfill(' ') << std::left << "[[ New memory high mark " << " : " << std::setw(8) << std::fixed << std::setprecision(2) << curr_mem_high_mark / (1024.0 * 1024.0) << " MB per task " << std::setw(23) << std::right << "]]" << std::endl;
global_mem_high_mark = curr_mem_high_mark;
}
}

View file

@ -97,7 +97,7 @@ int run( config_file& the_config )
//--------------------------------------------------------------------------------------------------------
//! number of resolution elements per dimension
const size_t ngrid = the_config.get_value<size_t>("setup", "GridRes");
music::COLUMN_ilog( "Grid resolution", ngrid );
//--------------------------------------------------------------------------------------------------------
//! box side length in h-1 Mpc
const real_t boxlen = the_config.get_value<double>("setup", "BoxLength");
@ -135,23 +135,28 @@ int run( config_file& the_config )
: ((lattice_str=="masked")? particle::lattice_masked
: particle::lattice_sc)))));
music::ilog << "Using " << lattice_str << " lattice for particle load." << std::endl;
// music::ilog << "Using " << lattice_str << " lattice for particle load." << std::endl;
music::COLUMN_ilog( "Lattice for particle load", lattice_str );
//--------------------------------------------------------------------------------------------------------
//! apply fixing of the complex mode amplitude following Angulo & Pontzen (2016) [https://arxiv.org/abs/1603.05253]
const bool bDoFixing = the_config.get_value_safe<bool>("setup", "DoFixing", false);
music::ilog << "Fixing of complex mode amplitudes is " << (bDoFixing?"enabled":"disabled") << std::endl;
music::COLUMN_ilog( "Fixing of mode amplitudes", bDoFixing );
// music::ilog << "Fixing of complex mode amplitudes is " << (bDoFixing?"enabled":"disabled") << std::endl;
const bool bDoInversion = the_config.get_value_safe<bool>("setup", "DoInversion", false);
music::ilog << "Inversion of the phase field is " << (bDoInversion?"enabled":"disabled") << std::endl;
music::COLUMN_ilog( "Inversion of the phase field", bDoInversion );
// music::ilog << "Inversion of the phase field is " << (bDoInversion?"enabled":"disabled") << std::endl;
//--------------------------------------------------------------------------------------------------------
//! do baryon ICs?
const bool bDoBaryons = the_config.get_value_safe<bool>("setup", "DoBaryons", false );
music::ilog << "Baryon ICs are " << (bDoBaryons?"enabled":"disabled") << std::endl;
// music::ilog << "Baryon ICs are " << (bDoBaryons?"enabled":"disabled") << std::endl;
music::COLUMN_ilog( "Baryon ICs", bDoBaryons );
//! enable also back-scaled decaying relative velocity mode? only first order!
const bool bDoLinearBCcorr = the_config.get_value_safe<bool>("setup", "DoBaryonVrel", false);
music::ilog << "Baryon linear relative velocity mode is " << (bDoLinearBCcorr?"enabled":"disabled") << std::endl;
// music::ilog << "Baryon linear relative velocity" << (bDoLinearBCcorr?"enabled":"disabled") << std::endl;
music::COLUMN_ilog( "Baryon linear relative velocity", bDoLinearBCcorr );
// compute mass fractions
std::map< cosmo_species, double > Omega;
if( bDoBaryons ){