14#include <deal.II/base/logstream.h>
15#include <deal.II/base/work_stream.h>
16#include <deal.II/numerics/vector_tools.h>
17#include <deal.II/numerics/vector_tools.templates.h>
24using namespace dealii;
28 template <
typename Description,
int dim,
typename Number>
30 : ParameterAcceptor(
"/A - TimeLoop")
31 , mpi_ensemble_(mpi_comm)
32 , hyperbolic_system_(mpi_ensemble_,
"/B - Equation")
33 , parabolic_system_(mpi_ensemble_,
"/B - Equation")
34 , discretization_(mpi_ensemble_,
"/C - Discretization")
35 , offline_data_(mpi_ensemble_, discretization_,
"/D - OfflineData")
36 , initial_values_(mpi_ensemble_,
42 , hyperbolic_module_(mpi_ensemble_,
47 "/F - HyperbolicModule")
48 , parabolic_module_(mpi_ensemble_,
54 "/G - ParabolicModule")
55 , time_integrator_(mpi_ensemble_,
59 "/H - TimeIntegrator")
60 , mesh_adaptor_(mpi_ensemble_,
64 hyperbolic_module_.initial_precomputed(),
65 hyperbolic_module_.alpha(),
67 , solution_transfer_(mpi_ensemble_,
72 , postprocessor_(mpi_ensemble_,
77 , vtu_output_(mpi_ensemble_,
82 hyperbolic_module_.initial_precomputed(),
83 hyperbolic_module_.alpha(),
86 , quantities_(mpi_ensemble_,
93 add_parameter(
"basename", base_name_,
"Base name for all output files");
95 t_final_ = Number(5.);
96 add_parameter(
"final time", t_final_,
"Final time");
98 enforce_t_final_ =
false;
99 add_parameter(
"enforce final time",
101 "Boolean indicating whether the final time should be "
102 "enforced strictly. If set to true the last time step is "
103 "shortened so that the simulation ends precisely at t_final");
105 timer_granularity_ = Number(0.01);
106 add_parameter(
"timer granularity",
108 "The timer granularity specifies the time interval after "
109 "which compute, output, postprocessing, and mesh adaptation "
110 "routines are run. This \"baseline tick\" is further "
111 "modified by the corresponding \"*_multiplier\" options");
113 enable_output_full_ =
false;
114 add_parameter(
"enable output full",
116 "Write out full pvtu records. The frequency is determined by "
117 "\"timer granularity\" and \"timer output full multiplier\"");
119 enable_output_levelsets_ =
false;
121 "enable output levelsets",
122 enable_output_levelsets_,
123 "Write out levelsets pvtu records. The frequency is determined by "
124 "\"timer granularity\" and \"timer output levelsets multiplier\"");
126 enable_compute_error_ =
false;
127 add_parameter(
"enable compute error",
128 enable_compute_error_,
129 "Flag to control whether we compute the Linfty Linf_norm of "
130 "the difference to an analytic solution. Implemented only "
131 "for certain initial state configurations.");
133 enable_compute_quantities_ =
false;
135 "enable compute quantities",
136 enable_compute_quantities_,
137 "Flag to control whether we compute quantities of interest. The "
138 "frequency how often quantities are logged is determined by \"timer "
139 "granularity\" and \"timer compute quantities multiplier\"");
141 enable_mesh_adaptivity_ =
false;
143 "enable mesh adaptivity",
144 enable_mesh_adaptivity_,
145 "Flag to control whether we use an adaptive mesh refinement strategy. "
146 "The frequency how often we query MeshAdaptor::analyze() for deciding "
147 "on adapting the mesh is determined by \"timer granularity\" and "
148 "\"timer mesh refinement multiplier\"");
150 timer_output_full_multiplier_ = 1;
151 add_parameter(
"timer output full multiplier",
152 timer_output_full_multiplier_,
153 "Multiplicative modifier applied to \"timer granularity\" "
154 "that determines the full pvtu writeout granularity");
156 timer_output_levelsets_multiplier_ = 1;
157 add_parameter(
"timer output levelsets multiplier",
158 timer_output_levelsets_multiplier_,
159 "Multiplicative modifier applied to \"timer granularity\" "
160 "that determines the levelsets pvtu writeout granularity");
162 timer_compute_quantities_multiplier_ = 1;
164 "timer compute quantities multiplier",
165 timer_compute_quantities_multiplier_,
166 "Multiplicative modifier applied to \"timer granularity\" that "
167 "determines the writeout granularity for quantities of interest");
169 std::copy(std::begin(View::component_names),
170 std::end(View::component_names),
171 std::back_inserter(error_quantities_));
173 add_parameter(
"error quantities",
175 "List of conserved quantities used in the computation of the "
178 error_normalize_ =
true;
179 add_parameter(
"error normalize",
181 "Flag to control whether the error should be normalized by "
182 "the corresponding norm of the analytic solution.");
185 add_parameter(
"resume", resume_,
"Resume an interrupted computation");
187 resume_at_time_zero_ =
false;
188 add_parameter(
"resume at time zero",
189 resume_at_time_zero_,
190 "Resume from the latest checkpoint but set the time to t=0.");
192 terminal_update_interval_ = 5;
193 add_parameter(
"terminal update interval",
194 terminal_update_interval_,
195 "Number of seconds after which output statistics are "
196 "recomputed and printed on the terminal. Setting the "
197 "interval to zero disables terminal output.");
199 terminal_correct_for_hypertreadhing_ =
true;
201 "terminal correct for hyperthreading",
202 terminal_correct_for_hypertreadhing_,
203 "If set to true, the CPU throughput is corrected by dividing the total "
204 "consumed CPU time by a factor of 2. This correction is only active if "
205 "the number of threads (per MPI rank) is 2.");
207 checkpoint_update_interval_ = 0;
209 "checkpoint update interval",
210 checkpoint_update_interval_,
211 "Number of seconds after which a new checkpoint is written out to "
212 "disk. Setting the interval to zero disables checkpointing.");
215 add_parameter(
"debug command",
217 "If set to a nonempty string then the host environment's "
218 "command processor is invoked via std::system() with the "
219 "specified string as command parameter.");
221 debug_filename_ =
"";
222 add_parameter(
"debug filename",
224 "If set to a nonempty string then we output the contents of "
225 "this file at the end. This is mainly useful in the "
226 "testsuite to output files we wish to compare");
237 template <
typename Description,
int dim,
typename Number>
241 std::cout <<
"TimeLoop<dim, Number>::run()" << std::endl;
245 base_name_ensemble_ = base_name_;
246 if (mpi_ensemble_.n_ensembles() > 1) {
247 print_info(
"setting up MPI ensemble");
248 base_name_ensemble_ +=
"-ensemble_" + dealii::Utilities::int_to_string(
249 mpi_ensemble_.ensemble(),
250 mpi_ensemble_.n_ensembles());
256 if (mpi_ensemble_.world_rank() == 0)
257 logfile_.open(base_name_ +
".log");
259 print_parameters(logfile_);
266 unsigned int timer_cycle = 0;
270 const auto prepare_compute_kernels = [&]() {
271 print_info(
"preparing compute kernels");
273 unsigned int n_parabolic_state_vectors =
274 parabolic_system_.get().n_parabolic_state_vectors();
276 offline_data_.prepare(
277 problem_dimension, n_precomputed_values, n_parabolic_state_vectors);
279 hyperbolic_module_.prepare();
280 parabolic_module_.prepare();
281 time_integrator_.prepare();
282 mesh_adaptor_.prepare( t);
283 postprocessor_.prepare();
284 vtu_output_.prepare();
285 quantities_.prepare(base_name_ensemble_);
286 print_mpi_partition(logfile_);
288 if (mpi_ensemble_.ensemble_rank() == 0)
289 n_global_dofs_ = dealii::Utilities::MPI::sum(
290 offline_data_.dof_handler().n_dofs(),
291 mpi_ensemble_.ensemble_leader_communicator());
295 Scope scope(computing_timer_,
"(re)initialize data structures");
296 print_info(
"initializing data structures");
299 print_info(
"resume: reading mesh and loading state vector");
301 read_checkpoint(state_vector,
305 prepare_compute_kernels);
307 if (resume_at_time_zero_) {
314 print_info(
"creating mesh and interpolating initial values");
316 discretization_.prepare(base_name_ensemble_);
318 prepare_compute_kernels();
320 Vectors::reinit_state_vector<Description>(state_vector, offline_data_);
321 std::get<0>(state_vector) =
322 initial_values_.get().interpolate_hyperbolic_vector();
324 Vectors::debug_poison_constrained_dofs<Description>(state_vector,
326 Vectors::debug_poison_precomputed_values<Description>(state_vector,
332 time_integrator_.prepare_state_vector(state_vector, t);
338 Number last_terminal_output = terminal_update_interval_ == Number(0.)
339 ? std::numeric_limits<Number>::max()
341 Number last_checkpoint = checkpoint_update_interval_ == Number(0.)
342 ? std::numeric_limits<Number>::max()
345 print_info(
"entering main loop");
346 computing_timer_[
"time loop"].start();
348 constexpr Number relax =
349 Number(1.) - Number(10.) * std::numeric_limits<Number>::epsilon();
351 unsigned int cycle = 1;
355 std::cout <<
"\n\n### cycle = " << cycle <<
" ###\n\n" << std::endl;
360 if (enable_compute_quantities_) {
361 Scope scope(computing_timer_,
362 "time step [X] - accumulate quantities");
363 quantities_.accumulate(state_vector, t);
368 if (t >= relax * timer_cycle * timer_granularity_) {
369 if (enable_compute_error_) {
378 Scope scope(computing_timer_,
379 "time step [X] - interpolate analytic solution");
380 Vectors::reinit_state_vector<Description>(analytic, offline_data_);
381 std::get<0>(analytic) =
382 initial_values_.get().interpolate_hyperbolic_vector(t);
386 base_name_ensemble_ +
"-analytic_solution",
391 output(state_vector, base_name_ensemble_ +
"-solution", t, timer_cycle);
393 if (enable_compute_quantities_ &&
394 (timer_cycle % timer_compute_quantities_multiplier_ == 0)) {
395 Scope scope(computing_timer_,
396 "time step [X] - write out quantities");
397 quantities_.write_out(state_vector, t, timer_cycle);
405 if (t >= relax * t_final_)
410 if (enable_mesh_adaptivity_) {
412 Scope scope(computing_timer_,
413 "time step [X] - analyze for mesh adaptation");
415 mesh_adaptor_.analyze(state_vector, t, cycle);
418 if (mesh_adaptor_.need_mesh_adaptation()) {
419 Scope scope_1(computing_timer_,
"(re)initialize data structures");
420 Scope scope_2(computing_timer_,
421 "time step [X] - perform mesh adaptation");
422 print_info(
"performing mesh adaptation");
424 adapt_mesh_and_transfer_state_vector(state_vector,
425 prepare_compute_kernels);
428 time_integrator_.prepare_state_vector(state_vector, t);
434 const auto tau = time_integrator_.step(
438 ? std::min(t_final_, timer_cycle * timer_granularity_)
439 : std::numeric_limits<Number>::max());
443 time_integrator_.prepare_state_vector(state_vector, t);
447 auto wall_time = computing_timer_[
"time loop"].wall_time();
449 Scope scope(computing_timer_,
450 "time step [X] _ - synchronization barriers");
452 Utilities::MPI::max(wall_time, mpi_ensemble_.world_communicator());
457 const bool write_to_log_file =
458 (terminal_update_interval_ != Number(0.)) &&
459 (t >= relax * timer_cycle * timer_granularity_);
461 const bool update_terminal =
462 (wall_time >= last_terminal_output + terminal_update_interval_);
464 if (write_to_log_file || update_terminal) {
465 Scope scope(computing_timer_,
466 "time step [X] _ - synchronization barriers");
467 print_cycle_statistics(cycle,
472 last_terminal_output = wall_time;
475 const bool update_checkpoint =
476 (wall_time >= last_checkpoint + checkpoint_update_interval_);
478 if (update_checkpoint) {
479 Scope scop(computing_timer_,
"time step [X] - perform checkpointing");
481 print_info(
"scheduling checkpointing");
482 write_checkpoint(state_vector, base_name_ensemble_, t, timer_cycle);
483 last_checkpoint = wall_time;
490 if (checkpoint_update_interval_ != Number(0.)) {
491 Scope scope(computing_timer_,
"time step [X] - perform checkpointing");
493 print_info(
"scheduling checkpointing");
494 write_checkpoint(state_vector, base_name_ensemble_, t, timer_cycle);
497 computing_timer_[
"time loop"].stop();
499 if (terminal_update_interval_ != Number(0.)) {
501 print_cycle_statistics(cycle,
509 if (enable_compute_error_) {
511 compute_error(state_vector, t);
518 if (mpi_ensemble_.world_rank() == 0) {
519 if (debug_command_ !=
"") {
520 auto result [[maybe_unused]] = std::system(debug_command_.c_str());
523 if (debug_filename_ !=
"") {
524 std::ifstream f(debug_filename_);
526 std::cout << f.rdbuf();
531 CALLGRIND_DUMP_STATS;
543 template <
typename Description,
int dim,
typename Number>
544 template <
typename Callable>
547 const std::string &base_name,
549 unsigned int &timer_cycle,
550 const Callable &prepare_compute_kernels)
553 std::cout <<
"TimeLoop<dim, Number>::read_checkpoint()" << std::endl;
560#if DEAL_II_VERSION_GTE(9, 6, 0)
561 discretization_.refinement() = 0;
562 discretization_.prepare(base_name);
563 discretization_.triangulation().load(base_name +
"-checkpoint.mesh");
567 dealii::ExcMessage(
"write_checkpoint() is not available with "
568 "deal.II versions prior to 9.6.0"));
571 prepare_compute_kernels();
577 std::string name = base_name +
"-checkpoint";
579 unsigned int transfer_handle;
580 if (mpi_ensemble_.ensemble_rank() == 0) {
581 std::string meta = name +
".metadata";
583 std::ifstream file(meta, std::ios::binary);
584 boost::archive::binary_iarchive ia(file);
585 ia >> t >> timer_cycle >> transfer_handle;
589 if constexpr (std::is_same_v<Number, double>)
591 &t, 1, MPI_DOUBLE, 0, mpi_ensemble_.ensemble_communicator());
594 MPI_Bcast(&t, 1, MPI_FLOAT, 0, mpi_ensemble_.ensemble_communicator());
595 AssertThrowMPI(ierr);
597 ierr = MPI_Bcast(&timer_cycle,
601 mpi_ensemble_.ensemble_communicator());
602 AssertThrowMPI(ierr);
604 ierr = MPI_Bcast(&transfer_handle,
608 mpi_ensemble_.ensemble_communicator());
609 AssertThrowMPI(ierr);
613 Vectors::reinit_state_vector<Description>(state_vector, offline_data_);
615 solution_transfer_.set_handle(transfer_handle);
616 solution_transfer_.project(state_vector);
617 solution_transfer_.reset_handle();
624 Vectors::debug_poison_constrained_dofs<Description>(state_vector,
627 Vectors::debug_poison_precomputed_values<Description>(state_vector,
629 time_integrator_.prepare_state_vector(state_vector, t);
633 template <
typename Description,
int dim,
typename Number>
636 const std::string &base_name,
638 const unsigned int &timer_cycle)
641 std::cout <<
"TimeLoop<dim, Number>::write_checkpoint()" << std::endl;
644 solution_transfer_.prepare_projection(state_vector);
645 const auto transfer_handle = solution_transfer_.get_handle();
646 solution_transfer_.reset_handle();
648 std::string name = base_name +
"-checkpoint";
650 if (mpi_ensemble_.ensemble_rank() == 0) {
651 for (
const std::string suffix :
652 {
".mesh",
".mesh_fixed.data",
".mesh.info",
".metadata"})
653 if (std::filesystem::exists(name + suffix))
654 std::filesystem::rename(name + suffix, name + suffix +
"~");
657#if DEAL_II_VERSION_GTE(9, 6, 0)
658 const auto &triangulation = discretization_.triangulation();
659 triangulation.save(name +
".mesh");
663 dealii::ExcMessage(
"write_checkpoint() is not available with "
664 "deal.II versions prior to 9.6.0"));
671 if (mpi_ensemble_.ensemble_rank() == 0) {
672 std::string meta = name +
".metadata";
673 std::ofstream file(meta, std::ios::binary | std::ios::trunc);
674 boost::archive::binary_oarchive oa(file);
675 oa << t << timer_cycle << transfer_handle;
678 const int ierr = MPI_Barrier(mpi_ensemble_.ensemble_communicator());
679 AssertThrowMPI(ierr);
683 template <
typename Description,
int dim,
typename Number>
684 template <
typename Callable>
686 StateVector &state_vector,
const Callable &prepare_compute_kernels)
689 std::cout <<
"TimeLoop<dim, Number>::adapt_mesh_and_transfer_state_vector()"
693 AssertThrow(mpi_ensemble_.n_ensembles() == 1, dealii::ExcNotImplemented());
699 auto &triangulation = discretization_.triangulation();
700 mesh_adaptor_.mark_cells_for_coarsening_and_refinement(triangulation);
702 triangulation.prepare_coarsening_and_refinement();
704 solution_transfer_.prepare_projection(state_vector);
708 triangulation.execute_coarsening_and_refinement();
709 prepare_compute_kernels();
711 Vectors::reinit_state_vector<Description>(state_vector, offline_data_);
712 solution_transfer_.project(state_vector);
713 solution_transfer_.reset_handle();
719 Vectors::debug_poison_constrained_dofs<Description>(state_vector,
721 Vectors::debug_poison_precomputed_values<Description>(state_vector,
726 template <
typename Description,
int dim,
typename Number>
732 std::cout <<
"TimeLoop<dim, Number>::compute_error()" << std::endl;
735 Vector<Number> difference_per_cell(
736 discretization_.triangulation().n_active_cells());
738 Number linf_norm = 0.;
742 const auto analytic_U =
743 initial_values_.get().interpolate_hyperbolic_vector(t);
744 const auto &U = std::get<0>(state_vector);
748 analytic_component.reinit(offline_data_.scalar_partitioner());
749 error_component.reinit(offline_data_.scalar_partitioner());
752 for (
const auto &entry : error_quantities_) {
753 const auto &names = View::component_names;
754 const auto pos = std::find(std::begin(names), std::end(names), entry);
755 if (pos == std::end(names)) {
758 dealii::ExcMessage(
"Unknown component name »" + entry +
"«"));
762 const auto index = std::distance(std::begin(names), pos);
764 analytic_U.extract_component(analytic_component, index);
768 Number linf_norm_analytic = 0.;
769 Number l1_norm_analytic = 0.;
770 Number l2_norm_analytic = 0.;
772 if (error_normalize_) {
774 Utilities::MPI::max(analytic_component.linfty_norm(),
775 mpi_ensemble_.ensemble_communicator());
777 VectorTools::integrate_difference(
778 discretization_.mapping(),
779 offline_data_.dof_handler(),
781 Functions::ZeroFunction<dim, Number>(),
783 discretization_.quadrature_high_order(),
784 VectorTools::L1_norm);
787 Utilities::MPI::sum(difference_per_cell.l1_norm(),
788 mpi_ensemble_.ensemble_communicator());
790 VectorTools::integrate_difference(
791 discretization_.mapping(),
792 offline_data_.dof_handler(),
794 Functions::ZeroFunction<dim, Number>(),
796 discretization_.quadrature_high_order(),
797 VectorTools::L2_norm);
799 l2_norm_analytic = Number(std::sqrt(
800 Utilities::MPI::sum(
std::pow(difference_per_cell.l2_norm(), 2),
801 mpi_ensemble_.ensemble_communicator())));
806 U.extract_component(error_component, index);
808 offline_data_.affine_constraints().distribute(error_component);
809 error_component.update_ghost_values();
810 error_component -= analytic_component;
812 const Number linf_norm_error = Utilities::MPI::max(
813 error_component.linfty_norm(), mpi_ensemble_.ensemble_communicator());
815 VectorTools::integrate_difference(discretization_.mapping(),
816 offline_data_.dof_handler(),
818 Functions::ZeroFunction<dim, Number>(),
820 discretization_.quadrature_high_order(),
821 VectorTools::L1_norm);
823 const Number l1_norm_error = Utilities::MPI::sum(
824 difference_per_cell.l1_norm(), mpi_ensemble_.ensemble_communicator());
826 VectorTools::integrate_difference(discretization_.mapping(),
827 offline_data_.dof_handler(),
829 Functions::ZeroFunction<dim, Number>(),
831 discretization_.quadrature_high_order(),
832 VectorTools::L2_norm);
834 const Number l2_norm_error = Number(std::sqrt(
835 Utilities::MPI::sum(
std::pow(difference_per_cell.l2_norm(), 2),
836 mpi_ensemble_.ensemble_communicator())));
838 if (error_normalize_) {
839 linf_norm += linf_norm_error / linf_norm_analytic;
840 l1_norm += l1_norm_error / l1_norm_analytic;
841 l2_norm += l2_norm_error / l2_norm_analytic;
843 linf_norm += linf_norm_error;
844 l1_norm += l1_norm_error;
845 l2_norm += l2_norm_error;
849 if (mpi_ensemble_.ensemble_rank() != 0)
857 if (mpi_ensemble_.n_ensembles() > 1) {
858 linf_norm = Utilities::MPI::sum(
859 linf_norm, mpi_ensemble_.ensemble_leader_communicator());
860 l1_norm = Utilities::MPI::sum(
861 l1_norm, mpi_ensemble_.ensemble_leader_communicator());
862 l2_norm = Utilities::MPI::sum(
863 l2_norm, mpi_ensemble_.ensemble_leader_communicator());
866 if (mpi_ensemble_.world_rank() != 0)
869 logfile_ << std::endl <<
"Computed errors:" << std::endl << std::endl;
870 logfile_ << std::setprecision(16);
872 std::string description =
873 error_normalize_ ?
"Normalized consolidated" :
"Consolidated";
875 logfile_ << description +
" Linf, L1, and L2 errors at final time \n";
876 logfile_ << std::setprecision(16);
877 logfile_ <<
"#dofs = " << n_global_dofs_ << std::endl;
878 logfile_ <<
"t = " << t << std::endl;
879 logfile_ <<
"Linf = " << linf_norm << std::endl;
880 logfile_ <<
"L1 = " << l1_norm << std::endl;
881 logfile_ <<
"L2 = " << l2_norm << std::endl;
883 std::cout << description +
" Linf, L1, and L2 errors at final time \n";
884 std::cout << std::setprecision(16);
885 std::cout <<
"#dofs = " << n_global_dofs_ << std::endl;
886 std::cout <<
"t = " << t << std::endl;
887 std::cout <<
"Linf = " << linf_norm << std::endl;
888 std::cout <<
"L1 = " << l1_norm << std::endl;
889 std::cout <<
"L2 = " << l2_norm << std::endl;
893 template <
typename Description,
int dim,
typename Number>
895 const std::string &name,
897 const unsigned int cycle)
900 std::cout <<
"TimeLoop<dim, Number>::output(t = " << t <<
")" << std::endl;
903 const bool do_full_output =
904 (cycle % timer_output_full_multiplier_ == 0) && enable_output_full_;
905 const bool do_levelsets =
906 (cycle % timer_output_levelsets_multiplier_ == 0) &&
907 enable_output_levelsets_;
910 if (!(do_full_output || do_levelsets))
915 Scope scope(computing_timer_,
"time step [X] - perform vtu output");
916 print_info(
"scheduling output");
918 postprocessor_.compute(state_vector);
925 postprocessor_.reset_bounds();
928 mesh_adaptor_.compute_smoothness_indicators(state_vector);
930 vtu_output_.schedule_output(
931 state_vector, name, t, cycle, do_full_output, do_levelsets);
942 template <
typename Description,
int dim,
typename Number>
946 if (mpi_ensemble_.world_rank() != 0)
955 stream << std::endl <<
"Run time parameters:" << std::endl << std::endl;
956 ParameterAcceptor::prm.print_parameters(
957 stream, ParameterHandler::OutputStyle::ShortPRM);
962 std::ofstream output(base_name_ +
"-parameters.prm");
963 ParameterAcceptor::prm.print_parameters(output, ParameterHandler::ShortPRM);
967 template <
typename Description,
int dim,
typename Number>
978 std::vector<double> values = {
979 (double)offline_data_.n_export_indices(),
980 (double)offline_data_.n_locally_internal(),
981 (double)offline_data_.n_locally_owned(),
982 (double)offline_data_.n_locally_relevant(),
983 (double)offline_data_.n_export_indices() /
984 (double)offline_data_.n_locally_relevant(),
985 (double)offline_data_.n_locally_internal() /
986 (double)offline_data_.n_locally_relevant(),
987 (double)offline_data_.n_locally_owned() /
988 (double)offline_data_.n_locally_relevant()};
992 Utilities::MPI::min_max_avg(values, mpi_ensemble_.world_communicator());
994 if (mpi_ensemble_.world_rank() != 0)
997 std::ostringstream output;
1000 dealii::Utilities::needed_digits(mpi_ensemble_.n_world_ranks());
1002 const auto print_snippet = [&output, n](
const std::string &name,
1003 const auto &values) {
1004 output << name <<
": ";
1006 output << std::setw(9) << (
unsigned int)values.min
1007 <<
" [p" << std::setw(n) << values.min_index <<
"] "
1008 << std::setw(9) << (
unsigned int)values.avg <<
" "
1009 << std::setw(9) << (
unsigned int)values.max
1010 <<
" [p" << std::setw(n) << values.max_index <<
"]";
1014 const auto print_percentages = [&output, n](
const auto &percentages) {
1015 output << std::endl <<
" ";
1016 output <<
" (" << std::setw(3) << std::setprecision(2)
1017 << percentages.min * 100 <<
"% )"
1018 <<
" [p" << std::setw(n) << percentages.min_index <<
"] "
1019 <<
" (" << std::setw(3) << std::setprecision(2)
1020 << percentages.avg * 100 <<
"% )"
1022 <<
" (" << std::setw(3) << std::setprecision(2)
1023 << percentages.max * 100 <<
"% )"
1024 <<
" [p" << std::setw(n) << percentages.max_index <<
"]";
1027 output << std::endl << std::endl <<
"Partition: ";
1028 print_snippet(
"exp", data[0]);
1029 print_percentages(data[4]);
1031 output << std::endl <<
" ";
1032 print_snippet(
"int", data[1]);
1033 print_percentages(data[5]);
1035 output << std::endl <<
" ";
1036 print_snippet(
"own", data[2]);
1037 print_percentages(data[6]);
1039 output << std::endl <<
" ";
1040 print_snippet(
"rel", data[3]);
1042 stream << output.str() << std::endl;
1046 template <
typename Description,
int dim,
typename Number>
1049 if (mpi_ensemble_.world_rank() != 0)
1052 std::cout <<
"[INFO] " << header << std::endl;
1056 template <
typename Description,
int dim,
typename Number>
1059 const std::string &secondary,
1060 std::ostream &stream)
1062 if (mpi_ensemble_.world_rank() != 0)
1065 const int header_size = header.size();
1066 const auto padded_header =
1067 std::string(std::max(0, 34 - header_size) / 2,
' ') + header +
1068 std::string(std::max(0, 35 - header_size) / 2,
' ');
1070 const int secondary_size = secondary.size();
1071 const auto padded_secondary =
1072 std::string(std::max(0, 34 - secondary_size) / 2,
' ') + secondary +
1073 std::string(std::max(0, 35 - secondary_size) / 2,
' ');
1077 stream <<
" ####################################################\n";
1078 stream <<
" #########" << padded_header <<
"#########\n";
1079 stream <<
" #########" << padded_secondary <<
"#########\n";
1080 stream <<
" ####################################################\n";
1081 stream << std::endl;
1086 template <
typename Description,
int dim,
typename Number>
1088 unsigned int timer_cycle,
1089 Number last_checkpoint,
1090 std::ostream &stream,
1093 static const std::string vectorization_name = [] {
1094 constexpr auto width = VectorizedArray<Number>::size();
1100 result = std::to_string(width * 8 *
sizeof(Number)) +
"bit packed ";
1102 if constexpr (std::is_same_v<Number, double>)
1103 return result +
"double";
1104 else if constexpr (std::is_same_v<Number, float>)
1105 return result +
"float";
1110 stream <<
"Information: (HYP) " << hyperbolic_system_.get().problem_name;
1111 if constexpr (!ParabolicSystem::is_identity) {
1112 stream <<
"\n (PAR) " << parabolic_system_.get().problem_name;
1114 stream <<
"\n [" << base_name_ <<
"] ";
1115 if (mpi_ensemble_.n_ensembles() > 1) {
1116 stream << mpi_ensemble_.n_ensembles() <<
" ensembles ";
1119 << n_global_dofs_ <<
" Qdofs on "
1120 << mpi_ensemble_.n_world_ranks() <<
" ranks / "
1122 << omp_get_max_threads() <<
" omp"
1124 <<
"[openmp disabled]"
1126#ifdef WITH_DEAL_II_THREADS
1127 <<
" (" << MultithreadInfo::n_threads() <<
" dealii)"
1129 <<
" threads <" << vectorization_name <<
">\n";
1131 stream <<
" Last output cycle "
1133 <<
" at t = " << timer_granularity_ * (timer_cycle - 1)
1136 if (enable_output_full_)
1138 if (enable_output_levelsets_)
1139 stream <<
"levelsets ";
1140 if (enable_compute_quantities_)
1141 stream <<
"quantities ";
1145 if (checkpoint_update_interval_ != Number(0.)) {
1146 const auto wall_time =
1147 Utilities::MPI::min_max_avg(computing_timer_[
"time loop"].wall_time(),
1148 mpi_ensemble_.world_communicator());
1151 stream <<
" Last checkpoint at FINAL TIME\n";
1153 stream <<
" Last checkpoint at wall time "
1154 << std::setprecision(2) << std::fixed << last_checkpoint
1155 <<
"s (" << std::setprecision(0)
1156 << std::max(0., wall_time.max - last_checkpoint)
1157 <<
"s ago, interval " << checkpoint_update_interval_ <<
"s)\n";
1163 template <
typename Description,
int dim,
typename Number>
1165 std::ostream &stream)
1167 Utilities::System::MemoryStats stats;
1168 Utilities::System::get_memory_stats(stats);
1170 Utilities::MPI::MinMaxAvg data = Utilities::MPI::min_max_avg(
1171 stats.VmRSS / 1024., mpi_ensemble_.world_communicator());
1173 if (mpi_ensemble_.world_rank() != 0)
1176 std::ostringstream output;
1179 dealii::Utilities::needed_digits(mpi_ensemble_.n_world_ranks());
1181 output <<
"\nMemory: [MiB]"
1182 << std::setw(8) << data.min
1183 <<
" [p" << std::setw(n) << data.min_index <<
"] "
1184 << std::setw(8) << data.avg <<
" "
1185 << std::setw(8) << data.max
1186 <<
" [p" << std::setw(n) << data.max_index <<
"]";
1188 stream << output.str() << std::endl;
1192 template <
typename Description,
int dim,
typename Number>
1195 std::vector<std::ostringstream> output(computing_timer_.size());
1197 const auto equalize = [&]() {
1199 std::max_element(output.begin(),
1201 [](
const auto &left,
const auto &right) {
1202 return left.str().length() < right.str().length();
1204 const auto length = ptr->str().length();
1205 for (
auto &it : output)
1206 it << std::string(length - it.str().length() + 1,
' ');
1209 const auto print_wall_time = [&](
auto &timer,
auto &stream) {
1210 const auto wall_time = Utilities::MPI::min_max_avg(
1211 timer.wall_time(), mpi_ensemble_.world_communicator());
1213 constexpr auto eps = std::numeric_limits<double>::epsilon();
1218 const auto skew_negative = std::max(
1219 100. * (wall_time.min - wall_time.avg) / wall_time.avg - eps, -99.9);
1220 const auto skew_positive = std::min(
1221 100. * (wall_time.max - wall_time.avg) / wall_time.avg + eps, 99.9);
1223 stream << std::setprecision(2) << std::fixed << std::setw(9)
1224 << wall_time.avg <<
"s [sk: " << std::setprecision(1)
1225 << std::setw(5) << std::fixed << skew_negative <<
"%/"
1226 << std::setw(4) << std::fixed << skew_positive <<
"%]";
1228 dealii::Utilities::needed_digits(mpi_ensemble_.n_world_ranks());
1229 stream <<
" [p" << std::setw(n) << wall_time.min_index <<
"/"
1230 << wall_time.max_index <<
"]";
1233 const auto cpu_time_statistics =
1234 Utilities::MPI::min_max_avg(computing_timer_[
"time loop"].cpu_time(),
1235 mpi_ensemble_.world_communicator());
1236 const double total_cpu_time = cpu_time_statistics.sum;
1238 const auto print_cpu_time =
1239 [&](
auto &timer,
auto &stream,
bool percentage) {
1240 const auto cpu_time = Utilities::MPI::min_max_avg(
1241 timer.cpu_time(), mpi_ensemble_.world_communicator());
1243 stream << std::setprecision(2) << std::fixed << std::setw(12)
1244 << cpu_time.sum <<
"s ";
1247 stream <<
"(" << std::setprecision(1) << std::setw(4)
1248 << 100. * cpu_time.sum / total_cpu_time <<
"%)";
1251 auto jt = output.begin();
1252 for (
auto &it : computing_timer_)
1253 *jt++ <<
" " << it.first;
1256 jt = output.begin();
1257 for (
auto &it : computing_timer_)
1258 print_wall_time(it.second, *jt++);
1261 jt = output.begin();
1262 bool compute_percentages =
false;
1263 for (
auto &it : computing_timer_) {
1264 print_cpu_time(it.second, *jt++, compute_percentages);
1265 if (it.first.starts_with(
"time loop"))
1266 compute_percentages =
true;
1270 if (mpi_ensemble_.world_rank() != 0)
1273 stream << std::endl <<
"Timer statistics:\n";
1274 for (
auto &it : output)
1275 stream << it.str() << std::endl;
1279 template <
typename Description,
int dim,
typename Number>
1281 unsigned int cycle, Number t, std::ostream &stream,
bool final_time)
1287 static struct Data {
1288 unsigned int cycle = 0;
1290 double cpu_time_sum = 0.;
1291 double cpu_time_avg = 0.;
1292 double cpu_time_min = 0.;
1293 double cpu_time_max = 0.;
1294 double wall_time = 0.;
1295 } previous, current;
1297 static double time_per_second_exp = 0.;
1304 current.cycle = cycle;
1307 const auto wall_time_statistics =
1308 Utilities::MPI::min_max_avg(computing_timer_[
"time loop"].wall_time(),
1309 mpi_ensemble_.world_communicator());
1310 current.wall_time = wall_time_statistics.max;
1312 const auto cpu_time_statistics =
1313 Utilities::MPI::min_max_avg(computing_timer_[
"time loop"].cpu_time(),
1314 mpi_ensemble_.world_communicator());
1315 current.cpu_time_sum = cpu_time_statistics.sum;
1316 current.cpu_time_avg = cpu_time_statistics.avg;
1317 current.cpu_time_min = cpu_time_statistics.min;
1318 current.cpu_time_max = cpu_time_statistics.max;
1326 double delta_cycles = current.cycle - previous.cycle;
1327 const double cycles_per_second =
1328 delta_cycles / (current.wall_time - previous.wall_time);
1330 const auto efficiency = time_integrator_.efficiency();
1331 const auto n_dofs =
static_cast<double>(n_global_dofs_);
1333 double wall_m_dofs_per_sec = delta_cycles * n_dofs * efficiency / 1.e6 /
1334 (current.wall_time - previous.wall_time);
1336 double cpu_m_dofs_per_sec = delta_cycles * n_dofs * efficiency / 1.e6 /
1337 (current.cpu_time_sum - previous.cpu_time_sum);
1340 const bool fudge_cpu_timings = terminal_correct_for_hypertreadhing_ &&
1342 (omp_get_max_threads() == 2);
1347 if (fudge_cpu_timings)
1348 cpu_m_dofs_per_sec *= 2.;
1350 double cpu_time_skew = (current.cpu_time_max - current.cpu_time_min -
1351 previous.cpu_time_max + previous.cpu_time_min) /
1354 cpu_time_skew = std::max(0., cpu_time_skew);
1356 const double cpu_time_skew_percentage =
1357 cpu_time_skew * delta_cycles /
1358 (current.cpu_time_avg - previous.cpu_time_avg);
1360 const double delta_time =
1361 (current.t - previous.t) / (current.cycle - previous.cycle);
1362 const double time_per_second =
1363 (current.t - previous.t) / (current.wall_time - previous.wall_time);
1367 std::ostringstream output;
1370 output << std::endl;
1372 output <<
"Throughput:\n "
1373 << (fudge_cpu_timings ?
"CPU*: " :
"CPU : ")
1374 << std::setprecision(4) << std::fixed << cpu_m_dofs_per_sec
1376 << std::scientific << 1. / cpu_m_dofs_per_sec * 1.e-6
1377 <<
" s/Qdof/substep)" << std::endl;
1379 output <<
" [cpu time skew: "
1380 << std::setprecision(2) << std::scientific << cpu_time_skew
1382 << std::setprecision(1) << std::setw(4) << std::setfill(
' ') << std::fixed
1383 << 100. * cpu_time_skew_percentage
1384 <<
"%)]" << std::endl;
1387 << std::setprecision(4) << std::fixed << wall_m_dofs_per_sec
1389 << std::scientific << 1. / wall_m_dofs_per_sec * 1.e-6
1390 <<
" s/Qdof/substep) ("
1391 << std::setprecision(2) << std::fixed << cycles_per_second
1392 <<
" cycles/s)" << std::endl;
1394 const auto &scheme = time_integrator_.time_stepping_scheme();
1396 << Patterns::Tools::Convert<TimeSteppingScheme>::to_string(scheme)
1398 << std::setprecision(2) << std::fixed << hyperbolic_module_.cfl()
1400 << std::setprecision(0) << std::fixed << hyperbolic_module_.n_restarts()
1402 << std::setprecision(0) << std::fixed << parabolic_module_.n_restarts()
1404 << std::setprecision(0) << std::fixed << hyperbolic_module_.n_warnings()
1406 << std::setprecision(0) << std::fixed << parabolic_module_.n_warnings()
1408 << std::setprecision(0) << std::fixed << hyperbolic_module_.n_corrections()
1410 << std::setprecision(0) << std::fixed << parabolic_module_.n_corrections()
1411 <<
" corr) ]" << std::endl;
1413 if constexpr (!ParabolicSystem::is_identity)
1414 parabolic_module_.print_solver_statistics(output);
1416 output <<
" [ dt = "
1417 << std::scientific << std::setprecision(2) << delta_time
1420 <<
" dt/s) ]" << std::endl;
1425 time_per_second_exp = 0.8 * time_per_second_exp + 0.2 * time_per_second;
1426 auto eta =
static_cast<unsigned int>(std::max(t_final_ - t, Number(0.)) /
1427 time_per_second_exp);
1429 output <<
"\n ETA : ";
1431 const unsigned int days = eta / (24 * 3600);
1433 output << days <<
" d ";
1437 const unsigned int hours = eta / 3600;
1439 output << hours <<
" h ";
1443 const unsigned int minutes = eta / 60;
1444 output << minutes <<
" min";
1446 output <<
" (terminal update every "
1447 << std::setprecision(2) << std::fixed << terminal_update_interval_
1450 if (mpi_ensemble_.world_rank() != 0)
1453 stream << output.str() << std::endl;
1457 template <
typename Description,
int dim,
typename Number>
1461 unsigned int timer_cycle,
1462 Number last_checkpoint,
1463 bool write_to_logfile,
1466 std::ostringstream output;
1470 std::ostringstream primary;
1472 primary <<
"FINAL (cycle " << Utilities::int_to_string(cycle, 6) <<
")";
1474 primary <<
"Cycle " << Utilities::int_to_string(cycle, 6)
1475 <<
" (" << std::fixed << std::setprecision(1)
1476 << t / t_final_ * 100 <<
"%)";
1479 std::ostringstream secondary;
1480 secondary <<
"at time t = " << std::setprecision(8) << std::fixed << t;
1482 print_head(primary.str(), secondary.str(), output);
1486 print_information(timer_cycle, last_checkpoint, output, final_time);
1487 print_memory_statistics(output);
1488 print_timers(output);
1489 print_throughput(cycle, t, output, final_time);
1492 if (mpi_ensemble_.world_rank() != 0)
1496 std::cout <<
"\033[2J\033[H";
1498 std::cout << output.str() << std::flush;
1500 if (write_to_logfile) {
1501 logfile_ <<
"\n" << output.str() << std::flush;
void write_checkpoint(const StateVector &state_vector, const std::string &base_name, const Number &t, const unsigned int &output_cycle)
Vectors::ScalarVector< Number > ScalarVector
void print_timers(std::ostream &stream)
void output(StateVector &state_vector, const std::string &name, const Number t, const unsigned int cycle)
void print_memory_statistics(std::ostream &stream)
void print_mpi_partition(std::ostream &stream)
void print_parameters(std::ostream &stream)
void compute_error(StateVector &state_vector, Number t)
void read_checkpoint(StateVector &state_vector, const std::string &base_name, Number &t, unsigned int &output_cycle, const Callable &prepare_compute_kernels)
TimeLoop(const MPI_Comm &mpi_comm)
void print_cycle_statistics(unsigned int cycle, Number t, unsigned int output_cycle, Number last_checkpoint, bool write_to_logfile=false, bool final_time=false)
typename View::StateVector StateVector
void print_head(const std::string &header, const std::string &secondary, std::ostream &stream)
void print_throughput(unsigned int cycle, Number t, std::ostream &stream, bool final_time=false)
void print_information(unsigned int output_cycle, Number last_checkpoint, std::ostream &stream, bool final_time=false)
void print_info(const std::string &header)
void adapt_mesh_and_transfer_state_vector(StateVector &state_vector, const Callable &prepare_compute_kernels)
T pow(const T x, const T b)
void print_revision_and_version(std::ostream &stream)