13#include <deal.II/base/logstream.h>
14#include <deal.II/base/work_stream.h>
15#include <deal.II/numerics/vector_tools.h>
16#include <deal.II/numerics/vector_tools.templates.h>
22using namespace dealii;
26 template <
typename Description,
int dim,
typename Number>
28 : ParameterAcceptor(
"/A - TimeLoop")
29 , mpi_ensemble_(mpi_comm)
30 , hyperbolic_system_(mpi_ensemble_,
"/B - Equation")
31 , parabolic_system_(mpi_ensemble_,
"/B - Equation")
32 , discretization_(mpi_ensemble_,
"/C - Discretization")
33 , offline_data_(mpi_ensemble_, discretization_,
"/D - OfflineData")
34 , initial_values_(mpi_ensemble_,
40 , hyperbolic_module_(mpi_ensemble_,
45 "/F - HyperbolicModule")
46 , parabolic_module_(mpi_ensemble_,
52 "/G - ParabolicModule")
53 , time_integrator_(mpi_ensemble_,
57 "/H - TimeIntegrator")
58 , mesh_adaptor_(mpi_ensemble_,
62 hyperbolic_module_.initial_precomputed(),
63 hyperbolic_module_.alpha(),
66 mpi_ensemble_, offline_data_, hyperbolic_system_, parabolic_system_)
67 , postprocessor_(mpi_ensemble_,
72 , vtu_output_(mpi_ensemble_,
77 hyperbolic_module_.initial_precomputed(),
78 hyperbolic_module_.alpha(),
80 , quantities_(mpi_ensemble_,
87 add_parameter(
"basename", base_name_,
"Base name for all output files");
89 t_final_ = Number(5.);
90 add_parameter(
"final time", t_final_,
"Final time");
92 enforce_t_final_ =
false;
93 add_parameter(
"enforce final time",
95 "Boolean indicating whether the final time should be "
96 "enforced strictly. If set to true the last time step is "
97 "shortened so that the simulation ends precisely at t_final");
99 timer_granularity_ = Number(0.01);
100 add_parameter(
"timer granularity",
102 "The timer granularity specifies the time interval after "
103 "which compute, output, postprocessing, and mesh adaptation "
104 "routines are run. This \"baseline tick\" is further "
105 "modified by the corresponding \"*_multiplier\" options");
107 enable_checkpointing_ =
false;
109 "enable checkpointing",
110 enable_checkpointing_,
111 "Write out checkpoints to resume an interrupted computation at timer "
112 "granularity intervals. The frequency is determined by \"timer "
113 "granularity\" and \"timer checkpoint multiplier\"");
115 enable_output_full_ =
false;
116 add_parameter(
"enable output full",
118 "Write out full pvtu records. The frequency is determined by "
119 "\"timer granularity\" and \"timer output full multiplier\"");
121 enable_output_levelsets_ =
false;
123 "enable output levelsets",
124 enable_output_levelsets_,
125 "Write out levelsets pvtu records. The frequency is determined by "
126 "\"timer granularity\" and \"timer output levelsets multiplier\"");
128 enable_compute_error_ =
false;
129 add_parameter(
"enable compute error",
130 enable_compute_error_,
131 "Flag to control whether we compute the Linfty Linf_norm of "
132 "the difference to an analytic solution. Implemented only "
133 "for certain initial state configurations.");
135 enable_compute_quantities_ =
false;
137 "enable compute quantities",
138 enable_compute_quantities_,
139 "Flag to control whether we compute quantities of interest. The "
140 "frequency how often quantities are logged is determined by \"timer "
141 "granularity\" and \"timer compute quantities multiplier\"");
143 enable_mesh_adaptivity_ =
false;
145 "enable mesh adaptivity",
146 enable_mesh_adaptivity_,
147 "Flag to control whether we use an adaptive mesh refinement strategy. "
148 "The frequency how often we query MeshAdaptor::analyze() for deciding "
149 "on adapting the mesh is determined by \"timer granularity\" and "
150 "\"timer mesh refinement multiplier\"");
152 timer_checkpoint_multiplier_ = 1;
153 add_parameter(
"timer checkpoint multiplier",
154 timer_checkpoint_multiplier_,
155 "Multiplicative modifier applied to \"timer granularity\" "
156 "that determines the checkpointing granularity");
158 timer_output_full_multiplier_ = 1;
159 add_parameter(
"timer output full multiplier",
160 timer_output_full_multiplier_,
161 "Multiplicative modifier applied to \"timer granularity\" "
162 "that determines the full pvtu writeout granularity");
164 timer_output_levelsets_multiplier_ = 1;
165 add_parameter(
"timer output levelsets multiplier",
166 timer_output_levelsets_multiplier_,
167 "Multiplicative modifier applied to \"timer granularity\" "
168 "that determines the levelsets pvtu writeout granularity");
170 timer_compute_quantities_multiplier_ = 1;
172 "timer compute quantities multiplier",
173 timer_compute_quantities_multiplier_,
174 "Multiplicative modifier applied to \"timer granularity\" that "
175 "determines the writeout granularity for quantities of interest");
177 std::copy(std::begin(View::component_names),
178 std::end(View::component_names),
179 std::back_inserter(error_quantities_));
181 add_parameter(
"error quantities",
183 "List of conserved quantities used in the computation of the "
186 error_normalize_ =
true;
187 add_parameter(
"error normalize",
189 "Flag to control whether the error should be normalized by "
190 "the corresponding norm of the analytic solution.");
193 add_parameter(
"resume", resume_,
"Resume an interrupted computation");
195 resume_at_time_zero_ =
false;
196 add_parameter(
"resume at time zero",
197 resume_at_time_zero_,
198 "Resume from the latest checkpoint but set the time to t=0.");
200 terminal_update_interval_ = 5;
201 add_parameter(
"terminal update interval",
202 terminal_update_interval_,
203 "Number of seconds after which output statistics are "
204 "recomputed and printed on the terminal");
206 terminal_show_rank_throughput_ =
true;
207 add_parameter(
"terminal show rank throughput",
208 terminal_show_rank_throughput_,
209 "If set to true an average per rank throughput is computed "
210 "by dividing the total consumed CPU time (per rank) by the "
211 "number of threads (per rank). If set to false then a plain "
212 "average per thread \"CPU\" throughput value is computed by "
213 "using the umodified total accumulated CPU time.");
215 debug_filename_ =
"";
216 add_parameter(
"debug filename",
218 "If set to a nonempty string then we output the contents of "
219 "this file at the end. This is mainly useful in the "
220 "testsuite to output files we wish to compare");
224 template <
typename Description,
int dim,
typename Number>
228 std::cout <<
"TimeLoop<dim, Number>::run()" << std::endl;
232 base_name_ensemble_ = base_name_;
233 if (mpi_ensemble_.n_ensembles() > 1) {
234 print_info(
"setting up MPI ensemble");
235 base_name_ensemble_ +=
"-ensemble_" + dealii::Utilities::int_to_string(
236 mpi_ensemble_.ensemble(),
237 mpi_ensemble_.n_ensembles());
243 if (mpi_ensemble_.world_rank() == 0)
244 logfile_.open(base_name_ +
".log");
246 print_parameters(logfile_);
253 unsigned int timer_cycle = 0;
257 const auto prepare_compute_kernels = [&]() {
258 print_info(
"preparing compute kernels");
260 unsigned int n_parabolic_state_vectors =
261 parabolic_system_.get().n_parabolic_state_vectors();
263 offline_data_.prepare(
264 problem_dimension, n_precomputed_values, n_parabolic_state_vectors);
266 hyperbolic_module_.prepare();
267 parabolic_module_.prepare();
268 time_integrator_.prepare();
269 mesh_adaptor_.prepare( t);
270 postprocessor_.prepare();
271 vtu_output_.prepare();
272 quantities_.prepare(base_name_ensemble_);
273 print_mpi_partition(logfile_);
275 if (mpi_ensemble_.ensemble_rank() == 0)
276 n_global_dofs_ = dealii::Utilities::MPI::sum(
277 offline_data_.dof_handler().n_dofs(),
278 mpi_ensemble_.ensemble_leader_communicator());
282 Scope scope(computing_timer_,
"(re)initialize data structures");
283 print_info(
"initializing data structures");
286 print_info(
"resume: reading mesh and loading state vector");
288 read_checkpoint(state_vector,
292 prepare_compute_kernels);
294 if (resume_at_time_zero_) {
301 print_info(
"creating mesh and interpolating initial values");
303 discretization_.prepare(base_name_ensemble_);
305 prepare_compute_kernels();
307 Vectors::reinit_state_vector<Description>(state_vector, offline_data_);
308 std::get<0>(state_vector) =
309 initial_values_.get().interpolate_hyperbolic_vector();
317 Vectors::debug_poison_constrained_dofs<Description>(state_vector,
319 Vectors::debug_poison_precomputed_values<Description>(state_vector,
322 unsigned int cycle = 1;
323 Number last_terminal_output = (terminal_update_interval_ == Number(0.)
324 ? std::numeric_limits<Number>::max()
325 : std::numeric_limits<Number>::lowest());
331 print_info(
"entering main loop");
332 computing_timer_[
"time loop"].start();
334 constexpr Number relax =
335 Number(1.) - Number(10.) * std::numeric_limits<Number>::epsilon();
340 std::cout <<
"\n\n### cycle = " << cycle <<
" ###\n\n" << std::endl;
345 if (enable_compute_quantities_) {
346 Scope scope(computing_timer_,
347 "time step [X] - accumulate quantities");
348 quantities_.accumulate(state_vector, t);
353 if (t >= relax * timer_cycle * timer_granularity_) {
354 if (enable_compute_error_) {
362 Scope scope(computing_timer_,
363 "time step [X] - interpolate analytic solution");
364 Vectors::reinit_state_vector<Description>(analytic, offline_data_);
365 std::get<0>(analytic) =
366 initial_values_.get().interpolate_hyperbolic_vector(t);
375 base_name_ensemble_ +
"-analytic_solution",
380 output(state_vector, base_name_ensemble_ +
"-solution", t, timer_cycle);
382 if (enable_compute_quantities_ &&
383 (timer_cycle % timer_compute_quantities_multiplier_ == 0)) {
384 Scope scope(computing_timer_,
385 "time step [X] - write out quantities");
386 quantities_.write_out(state_vector, t, timer_cycle);
394 if (t >= relax * t_final_)
399 if (enable_mesh_adaptivity_) {
401 Scope scope(computing_timer_,
402 "time step [X] - analyze for mesh adaptation");
404 mesh_adaptor_.analyze(state_vector, t, cycle);
407 if (mesh_adaptor_.need_mesh_adaptation()) {
408 Scope scope(computing_timer_,
"(re)initialize data structures");
409 print_info(
"performing mesh adaptation");
411 hyperbolic_module_.prepare_state_vector(state_vector, t);
412 adapt_mesh_and_transfer_state_vector(state_vector,
413 prepare_compute_kernels);
419 const auto tau = time_integrator_.step(
423 ? std::min(t_final_, timer_cycle * timer_granularity_)
424 : std::numeric_limits<Number>::max());
429 if (terminal_update_interval_ != Number(0.)) {
432 const bool write_to_log_file =
433 (t >= relax * timer_cycle * timer_granularity_);
436 const auto wall_time = computing_timer_[
"time loop"].wall_time();
437 int update_terminal =
438 (wall_time >= last_terminal_output + terminal_update_interval_);
441 const auto ierr = MPI_Bcast(&update_terminal,
445 mpi_ensemble_.world_communicator());
446 AssertThrowMPI(ierr);
448 if (write_to_log_file || update_terminal) {
449 print_cycle_statistics(
450 cycle, t, timer_cycle, write_to_log_file);
451 last_terminal_output = wall_time;
459 computing_timer_[
"time loop"].stop();
461 if (terminal_update_interval_ != Number(0.)) {
463 print_cycle_statistics(
464 cycle, t, timer_cycle,
true,
true);
467 if (enable_compute_error_) {
469 compute_error(state_vector, t);
472 if (mpi_ensemble_.world_rank() == 0 && debug_filename_ !=
"") {
473 std::ifstream f(debug_filename_);
475 std::cout << f.rdbuf();
479 CALLGRIND_DUMP_STATS;
484 template <
typename Description,
int dim,
typename Number>
485 template <
typename Callable>
488 const std::string &base_name,
490 unsigned int &output_cycle,
491 const Callable &prepare_compute_kernels)
494 std::cout <<
"TimeLoop<dim, Number>::read_checkpoint()" << std::endl;
497 AssertThrow(have_distributed_triangulation<dim>,
499 "read_checkpoint() is not implemented for "
500 "distributed::shared::Triangulation which we use in 1D"));
506#if !DEAL_II_VERSION_GTE(9, 6, 0)
507 if constexpr (have_distributed_triangulation<dim>) {
509 discretization_.refinement() = 0;
510 discretization_.prepare(base_name);
511 discretization_.triangulation().load(base_name +
"-checkpoint.mesh");
512#if !DEAL_II_VERSION_GTE(9, 6, 0)
516 prepare_compute_kernels();
522 std::string name = base_name +
"-checkpoint";
524 unsigned int transfer_handle;
525 if (mpi_ensemble_.ensemble_rank() == 0) {
526 std::string meta = name +
".metadata";
528 std::ifstream file(meta, std::ios::binary);
529 boost::archive::binary_iarchive ia(file);
530 ia >> t >> output_cycle >> transfer_handle;
534 if constexpr (std::is_same_v<Number, double>)
536 &t, 1, MPI_DOUBLE, 0, mpi_ensemble_.ensemble_communicator());
539 MPI_Bcast(&t, 1, MPI_FLOAT, 0, mpi_ensemble_.ensemble_communicator());
540 AssertThrowMPI(ierr);
542 ierr = MPI_Bcast(&output_cycle,
546 mpi_ensemble_.ensemble_communicator());
547 AssertThrowMPI(ierr);
549 ierr = MPI_Bcast(&transfer_handle,
553 mpi_ensemble_.ensemble_communicator());
554 AssertThrowMPI(ierr);
558 Vectors::reinit_state_vector<Description>(state_vector, offline_data_);
560 solution_transfer_.set_handle(transfer_handle);
561 solution_transfer_.project(state_vector);
562 solution_transfer_.reset_handle();
566 template <
typename Description,
int dim,
typename Number>
569 const std::string &base_name,
571 const unsigned int &output_cycle)
574 std::cout <<
"TimeLoop<dim, Number>::write_checkpoint()" << std::endl;
577 AssertThrow(have_distributed_triangulation<dim>,
579 "write_checkpoint() is not implemented for "
580 "distributed::shared::Triangulation which we use in 1D"));
583 solution_transfer_.prepare_projection(state_vector);
584 const auto transfer_handle = solution_transfer_.get_handle();
585 solution_transfer_.reset_handle();
587 std::string name = base_name +
"-checkpoint";
589 if (mpi_ensemble_.ensemble_rank() == 0) {
590 for (
const std::string suffix :
591 {
".mesh",
".mesh_fixed.data",
".mesh.info",
".metadata"})
592 if (std::filesystem::exists(name + suffix))
593 std::filesystem::rename(name + suffix, name + suffix +
"~");
596#if !DEAL_II_VERSION_GTE(9, 6, 0)
597 if constexpr (have_distributed_triangulation<dim>) {
599 const auto &triangulation = discretization_.triangulation();
600 triangulation.save(name +
".mesh");
601#if !DEAL_II_VERSION_GTE(9, 6, 0)
609 if (mpi_ensemble_.ensemble_rank() == 0) {
610 std::string meta = name +
".metadata";
611 std::ofstream file(meta, std::ios::binary | std::ios::trunc);
612 boost::archive::binary_oarchive oa(file);
613 oa << t << output_cycle << transfer_handle;
616 const int ierr = MPI_Barrier(mpi_ensemble_.ensemble_communicator());
617 AssertThrowMPI(ierr);
621 template <
typename Description,
int dim,
typename Number>
622 template <
typename Callable>
624 StateVector &state_vector,
const Callable &prepare_compute_kernels)
627 std::cout <<
"TimeLoop<dim, Number>::adapt_mesh_and_transfer_state_vector()"
631 AssertThrow(mpi_ensemble_.n_ensembles() == 1, dealii::ExcNotImplemented());
637 auto &triangulation = discretization_.triangulation();
638 mesh_adaptor_.mark_cells_for_coarsening_and_refinement(triangulation);
640 triangulation.prepare_coarsening_and_refinement();
643 solution_transfer_.prepare_projection(state_vector);
647 triangulation.execute_coarsening_and_refinement();
648 prepare_compute_kernels();
650 Vectors::reinit_state_vector<Description>(state_vector, offline_data_);
651 solution_transfer_.project(state_vector);
652 solution_transfer_.reset_handle();
656 template <
typename Description,
int dim,
typename Number>
662 std::cout <<
"TimeLoop<dim, Number>::compute_error()" << std::endl;
665 hyperbolic_module_.prepare_state_vector(state_vector, t);
667 Vector<Number> difference_per_cell(
668 discretization_.triangulation().n_active_cells());
670 Number linf_norm = 0.;
674 const auto analytic_U =
675 initial_values_.get().interpolate_hyperbolic_vector(t);
676 const auto &U = std::get<0>(state_vector);
680 analytic_component.reinit(offline_data_.scalar_partitioner());
681 error_component.reinit(offline_data_.scalar_partitioner());
684 for (
const auto &entry : error_quantities_) {
685 const auto &names = View::component_names;
686 const auto pos = std::find(std::begin(names), std::end(names), entry);
687 if (pos == std::end(names)) {
690 dealii::ExcMessage(
"Unknown component name »" + entry +
"«"));
694 const auto index = std::distance(std::begin(names), pos);
696 analytic_U.extract_component(analytic_component, index);
700 Number linf_norm_analytic = 0.;
701 Number l1_norm_analytic = 0.;
702 Number l2_norm_analytic = 0.;
704 if (error_normalize_) {
706 Utilities::MPI::max(analytic_component.linfty_norm(),
707 mpi_ensemble_.ensemble_communicator());
709 VectorTools::integrate_difference(
710 offline_data_.dof_handler(),
712 Functions::ZeroFunction<dim, Number>(),
715 VectorTools::L1_norm);
718 Utilities::MPI::sum(difference_per_cell.l1_norm(),
719 mpi_ensemble_.ensemble_communicator());
721 VectorTools::integrate_difference(
722 offline_data_.dof_handler(),
724 Functions::ZeroFunction<dim, Number>(),
727 VectorTools::L2_norm);
729 l2_norm_analytic = Number(std::sqrt(
730 Utilities::MPI::sum(
std::pow(difference_per_cell.l2_norm(), 2),
731 mpi_ensemble_.ensemble_communicator())));
736 U.extract_component(error_component, index);
738 offline_data_.affine_constraints().distribute(error_component);
739 error_component.update_ghost_values();
740 error_component -= analytic_component;
742 const Number linf_norm_error = Utilities::MPI::max(
743 error_component.linfty_norm(), mpi_ensemble_.ensemble_communicator());
745 VectorTools::integrate_difference(offline_data_.dof_handler(),
747 Functions::ZeroFunction<dim, Number>(),
750 VectorTools::L1_norm);
752 const Number l1_norm_error = Utilities::MPI::sum(
753 difference_per_cell.l1_norm(), mpi_ensemble_.ensemble_communicator());
755 VectorTools::integrate_difference(offline_data_.dof_handler(),
757 Functions::ZeroFunction<dim, Number>(),
760 VectorTools::L2_norm);
762 const Number l2_norm_error = Number(std::sqrt(
763 Utilities::MPI::sum(
std::pow(difference_per_cell.l2_norm(), 2),
764 mpi_ensemble_.ensemble_communicator())));
766 if (error_normalize_) {
767 linf_norm += linf_norm_error / linf_norm_analytic;
768 l1_norm += l1_norm_error / l1_norm_analytic;
769 l2_norm += l2_norm_error / l2_norm_analytic;
771 linf_norm += linf_norm_error;
772 l1_norm += l1_norm_error;
773 l2_norm += l2_norm_error;
777 if (mpi_ensemble_.ensemble_rank() != 0)
785 if (mpi_ensemble_.n_ensembles() > 1) {
786 linf_norm = Utilities::MPI::sum(
787 linf_norm, mpi_ensemble_.ensemble_leader_communicator());
788 l1_norm = Utilities::MPI::sum(
789 l1_norm, mpi_ensemble_.ensemble_leader_communicator());
790 l2_norm = Utilities::MPI::sum(
791 l2_norm, mpi_ensemble_.ensemble_leader_communicator());
794 if (mpi_ensemble_.world_rank() != 0)
797 logfile_ << std::endl <<
"Computed errors:" << std::endl << std::endl;
798 logfile_ << std::setprecision(16);
800 std::string description =
801 error_normalize_ ?
"Normalized consolidated" :
"Consolidated";
803 logfile_ << description +
" Linf, L1, and L2 errors at final time \n";
804 logfile_ << std::setprecision(16);
805 logfile_ <<
"#dofs = " << n_global_dofs_ << std::endl;
806 logfile_ <<
"t = " << t << std::endl;
807 logfile_ <<
"Linf = " << linf_norm << std::endl;
808 logfile_ <<
"L1 = " << l1_norm << std::endl;
809 logfile_ <<
"L2 = " << l2_norm << std::endl;
811 std::cout << description +
" Linf, L1, and L2 errors at final time \n";
812 std::cout << std::setprecision(16);
813 std::cout <<
"#dofs = " << n_global_dofs_ << std::endl;
814 std::cout <<
"t = " << t << std::endl;
815 std::cout <<
"Linf = " << linf_norm << std::endl;
816 std::cout <<
"L1 = " << l1_norm << std::endl;
817 std::cout <<
"L2 = " << l2_norm << std::endl;
821 template <
typename Description,
int dim,
typename Number>
823 const std::string &name,
825 const unsigned int cycle)
828 std::cout <<
"TimeLoop<dim, Number>::output(t = " << t <<
")" << std::endl;
831 const bool do_full_output =
832 (cycle % timer_output_full_multiplier_ == 0) && enable_output_full_;
833 const bool do_levelsets =
834 (cycle % timer_output_levelsets_multiplier_ == 0) &&
835 enable_output_levelsets_;
836 const bool do_checkpointing =
837 (cycle % timer_checkpoint_multiplier_ == 0) && enable_checkpointing_;
840 if (!(do_full_output || do_levelsets || do_checkpointing))
843 hyperbolic_module_.prepare_state_vector(state_vector, t);
846 if (do_full_output || do_levelsets) {
847 Scope scope(computing_timer_,
"time step [X] - perform vtu output");
848 print_info(
"scheduling output");
850 postprocessor_.compute(state_vector);
857 postprocessor_.reset_bounds();
859 vtu_output_.schedule_output(
860 state_vector, name, t, cycle, do_full_output, do_levelsets);
864 if (do_checkpointing) {
865 Scope scope(computing_timer_,
"time step [X] - perform checkpointing");
866 print_info(
"scheduling checkpointing");
867 write_checkpoint(state_vector, base_name_ensemble_, t, cycle);
877 template <
typename Description,
int dim,
typename Number>
881 if (mpi_ensemble_.world_rank() != 0)
890 stream << std::endl <<
"Run time parameters:" << std::endl << std::endl;
891 ParameterAcceptor::prm.print_parameters(
892 stream, ParameterHandler::OutputStyle::ShortPRM);
897 std::ofstream output(base_name_ +
"-parameters.prm");
898 ParameterAcceptor::prm.print_parameters(output, ParameterHandler::ShortPRM);
902 template <
typename Description,
int dim,
typename Number>
913 std::vector<double> values = {
914 (double)offline_data_.n_export_indices(),
915 (double)offline_data_.n_locally_internal(),
916 (double)offline_data_.n_locally_owned(),
917 (double)offline_data_.n_locally_relevant(),
918 (double)offline_data_.n_export_indices() /
919 (double)offline_data_.n_locally_relevant(),
920 (double)offline_data_.n_locally_internal() /
921 (double)offline_data_.n_locally_relevant(),
922 (double)offline_data_.n_locally_owned() /
923 (double)offline_data_.n_locally_relevant()};
927 Utilities::MPI::min_max_avg(values, mpi_ensemble_.world_communicator());
929 if (mpi_ensemble_.world_rank() != 0)
932 std::ostringstream output;
935 dealii::Utilities::needed_digits(mpi_ensemble_.n_world_ranks());
937 const auto print_snippet = [&output, n](
const std::string &name,
938 const auto &values) {
939 output << name <<
": ";
941 output << std::setw(9) << (
unsigned int)values.min
942 <<
" [p" << std::setw(n) << values.min_index <<
"] "
943 << std::setw(9) << (
unsigned int)values.avg <<
" "
944 << std::setw(9) << (
unsigned int)values.max
945 <<
" [p" << std::setw(n) << values.max_index <<
"]";
949 const auto print_percentages = [&output, n](
const auto &percentages) {
950 output << std::endl <<
" ";
951 output <<
" (" << std::setw(3) << std::setprecision(2)
952 << percentages.min * 100 <<
"% )"
953 <<
" [p" << std::setw(n) << percentages.min_index <<
"] "
954 <<
" (" << std::setw(3) << std::setprecision(2)
955 << percentages.avg * 100 <<
"% )"
957 <<
" (" << std::setw(3) << std::setprecision(2)
958 << percentages.max * 100 <<
"% )"
959 <<
" [p" << std::setw(n) << percentages.max_index <<
"]";
962 output << std::endl << std::endl <<
"Partition: ";
963 print_snippet(
"exp", data[0]);
964 print_percentages(data[4]);
966 output << std::endl <<
" ";
967 print_snippet(
"int", data[1]);
968 print_percentages(data[5]);
970 output << std::endl <<
" ";
971 print_snippet(
"own", data[2]);
972 print_percentages(data[6]);
974 output << std::endl <<
" ";
975 print_snippet(
"rel", data[3]);
977 stream << output.str() << std::endl;
981 template <
typename Description,
int dim,
typename Number>
983 std::ostream &stream)
985 Utilities::System::MemoryStats stats;
986 Utilities::System::get_memory_stats(stats);
988 Utilities::MPI::MinMaxAvg data = Utilities::MPI::min_max_avg(
989 stats.VmRSS / 1024., mpi_ensemble_.world_communicator());
991 if (mpi_ensemble_.world_rank() != 0)
994 std::ostringstream output;
997 dealii::Utilities::needed_digits(mpi_ensemble_.n_world_ranks());
999 output <<
"\nMemory: [MiB]"
1000 << std::setw(8) << data.min
1001 <<
" [p" << std::setw(n) << data.min_index <<
"] "
1002 << std::setw(8) << data.avg <<
" "
1003 << std::setw(8) << data.max
1004 <<
" [p" << std::setw(n) << data.max_index <<
"]";
1006 stream << output.str() << std::endl;
1010 template <
typename Description,
int dim,
typename Number>
1013 std::vector<std::ostringstream> output(computing_timer_.size());
1015 const auto equalize = [&]() {
1017 std::max_element(output.begin(),
1019 [](
const auto &left,
const auto &right) {
1020 return left.str().length() < right.str().length();
1022 const auto length = ptr->str().length();
1023 for (
auto &it : output)
1024 it << std::string(length - it.str().length() + 1,
' ');
1027 const auto print_wall_time = [&](
auto &timer,
auto &stream) {
1028 const auto wall_time = Utilities::MPI::min_max_avg(
1029 timer.wall_time(), mpi_ensemble_.world_communicator());
1031 constexpr auto eps = std::numeric_limits<double>::epsilon();
1036 const auto skew_negative = std::max(
1037 100. * (wall_time.min - wall_time.avg) / wall_time.avg - eps, -99.9);
1038 const auto skew_positive = std::min(
1039 100. * (wall_time.max - wall_time.avg) / wall_time.avg + eps, 99.9);
1041 stream << std::setprecision(2) << std::fixed << std::setw(8)
1042 << wall_time.avg <<
"s [sk: " << std::setprecision(1)
1043 << std::setw(5) << std::fixed << skew_negative <<
"%/"
1044 << std::setw(4) << std::fixed << skew_positive <<
"%]";
1046 dealii::Utilities::needed_digits(mpi_ensemble_.n_world_ranks());
1047 stream <<
" [p" << std::setw(n) << wall_time.min_index <<
"/"
1048 << wall_time.max_index <<
"]";
1051 const auto cpu_time_statistics =
1052 Utilities::MPI::min_max_avg(computing_timer_[
"time loop"].cpu_time(),
1053 mpi_ensemble_.world_communicator());
1054 const double total_cpu_time = cpu_time_statistics.sum;
1056 const auto print_cpu_time =
1057 [&](
auto &timer,
auto &stream,
bool percentage) {
1058 const auto cpu_time = Utilities::MPI::min_max_avg(
1059 timer.cpu_time(), mpi_ensemble_.world_communicator());
1061 stream << std::setprecision(2) << std::fixed << std::setw(9)
1062 << cpu_time.sum <<
"s ";
1065 stream <<
"(" << std::setprecision(1) << std::setw(4)
1066 << 100. * cpu_time.sum / total_cpu_time <<
"%)";
1069 auto jt = output.begin();
1070 for (
auto &it : computing_timer_)
1071 *jt++ <<
" " << it.first;
1074 jt = output.begin();
1075 for (
auto &it : computing_timer_)
1076 print_wall_time(it.second, *jt++);
1079 jt = output.begin();
1080 bool compute_percentages =
false;
1081 for (
auto &it : computing_timer_) {
1082 print_cpu_time(it.second, *jt++, compute_percentages);
1083 if (it.first.starts_with(
"time loop"))
1084 compute_percentages =
true;
1088 if (mpi_ensemble_.world_rank() != 0)
1091 stream << std::endl <<
"Timer statistics:\n";
1092 for (
auto &it : output)
1093 stream << it.str() << std::endl;
1097 template <
typename Description,
int dim,
typename Number>
1099 unsigned int cycle, Number t, std::ostream &stream,
bool final_time)
1105 static struct Data {
1106 unsigned int cycle = 0;
1108 double cpu_time_sum = 0.;
1109 double cpu_time_avg = 0.;
1110 double cpu_time_min = 0.;
1111 double cpu_time_max = 0.;
1112 double wall_time = 0.;
1113 } previous, current;
1115 static double time_per_second_exp = 0.;
1122 current.cycle = cycle;
1125 const auto wall_time_statistics =
1126 Utilities::MPI::min_max_avg(computing_timer_[
"time loop"].wall_time(),
1127 mpi_ensemble_.world_communicator());
1128 current.wall_time = wall_time_statistics.max;
1130 const auto cpu_time_statistics =
1131 Utilities::MPI::min_max_avg(computing_timer_[
"time loop"].cpu_time(),
1132 mpi_ensemble_.world_communicator());
1133 current.cpu_time_sum = cpu_time_statistics.sum;
1134 current.cpu_time_avg = cpu_time_statistics.avg;
1135 current.cpu_time_min = cpu_time_statistics.min;
1136 current.cpu_time_max = cpu_time_statistics.max;
1144 double delta_cycles = current.cycle - previous.cycle;
1145 const double cycles_per_second =
1146 delta_cycles / (current.wall_time - previous.wall_time);
1148 const auto efficiency = time_integrator_.efficiency();
1149 const auto n_dofs =
static_cast<double>(n_global_dofs_);
1151 const double wall_m_dofs_per_sec =
1152 delta_cycles * n_dofs / 1.e6 /
1153 (current.wall_time - previous.wall_time) * efficiency;
1155 double cpu_m_dofs_per_sec = delta_cycles * n_dofs / 1.e6 /
1156 (current.cpu_time_sum - previous.cpu_time_sum) *
1159 if (terminal_show_rank_throughput_)
1160 cpu_m_dofs_per_sec *= MultithreadInfo::n_threads();
1163 double cpu_time_skew = (current.cpu_time_max - current.cpu_time_min -
1164 previous.cpu_time_max + previous.cpu_time_min) /
1167 cpu_time_skew = std::max(0., cpu_time_skew);
1169 const double cpu_time_skew_percentage =
1170 cpu_time_skew * delta_cycles /
1171 (current.cpu_time_avg - previous.cpu_time_avg);
1173 const double delta_time =
1174 (current.t - previous.t) / (current.cycle - previous.cycle);
1175 const double time_per_second =
1176 (current.t - previous.t) / (current.wall_time - previous.wall_time);
1180 std::ostringstream output;
1183 output << std::endl;
1185 output <<
"Throughput:\n "
1186 << (terminal_show_rank_throughput_?
"RANK: " :
"CPU : ")
1187 << std::setprecision(4) << std::fixed << cpu_m_dofs_per_sec
1189 << std::scientific << 1. / cpu_m_dofs_per_sec * 1.e-6
1190 <<
" s/Qdof/substep)" << std::endl;
1192 output <<
" [cpu time skew: "
1193 << std::setprecision(2) << std::scientific << cpu_time_skew
1195 << std::setprecision(1) << std::setw(4) << std::setfill(
' ') << std::fixed
1196 << 100. * cpu_time_skew_percentage
1197 <<
"%)]" << std::endl;
1200 << std::setprecision(4) << std::fixed << wall_m_dofs_per_sec
1202 << std::scientific << 1. / wall_m_dofs_per_sec * 1.e-6
1203 <<
" s/Qdof/substep) ("
1204 << std::setprecision(2) << std::fixed << cycles_per_second
1205 <<
" cycles/s)" << std::endl;
1207 const auto &scheme = time_integrator_.time_stepping_scheme();
1209 << Patterns::Tools::Convert<TimeSteppingScheme>::to_string(scheme)
1211 << std::setprecision(2) << std::fixed << hyperbolic_module_.cfl()
1213 << std::setprecision(0) << std::fixed << hyperbolic_module_.n_restarts()
1215 << std::setprecision(0) << std::fixed << parabolic_module_.n_restarts()
1217 << std::setprecision(0) << std::fixed << hyperbolic_module_.n_warnings()
1219 << std::setprecision(0) << std::fixed << parabolic_module_.n_warnings()
1221 << std::setprecision(0) << std::fixed << hyperbolic_module_.n_corrections()
1223 << std::setprecision(0) << std::fixed << parabolic_module_.n_corrections()
1224 <<
" corr) ]" << std::endl;
1226 if constexpr (!ParabolicSystem::is_identity)
1227 parabolic_module_.print_solver_statistics(output);
1229 output <<
" [ dt = "
1230 << std::scientific << std::setprecision(2) << delta_time
1233 <<
" dt/s) ]" << std::endl;
1237 time_per_second_exp = 0.8 * time_per_second_exp + 0.2 * time_per_second;
1238 auto eta =
static_cast<unsigned int>(std::max(t_final_ - t, Number(0.)) /
1239 time_per_second_exp);
1241 output <<
"\n ETA : ";
1243 const unsigned int days = eta / (24 * 3600);
1245 output << days <<
" d ";
1249 const unsigned int hours = eta / 3600;
1251 output << hours <<
" h ";
1255 const unsigned int minutes = eta / 60;
1256 output << minutes <<
" min";
1258 if (mpi_ensemble_.world_rank() != 0)
1261 stream << output.str() << std::endl;
1265 template <
typename Description,
int dim,
typename Number>
1268 if (mpi_ensemble_.world_rank() != 0)
1271 std::cout <<
"[INFO] " << header << std::endl;
1275 template <
typename Description,
int dim,
typename Number>
1278 const std::string &secondary,
1279 std::ostream &stream)
1281 if (mpi_ensemble_.world_rank() != 0)
1284 const int header_size = header.size();
1285 const auto padded_header =
1286 std::string(std::max(0, 34 - header_size) / 2,
' ') + header +
1287 std::string(std::max(0, 35 - header_size) / 2,
' ');
1289 const int secondary_size = secondary.size();
1290 const auto padded_secondary =
1291 std::string(std::max(0, 34 - secondary_size) / 2,
' ') + secondary +
1292 std::string(std::max(0, 35 - secondary_size) / 2,
' ');
1296 stream <<
" ####################################################\n";
1297 stream <<
" #########" << padded_header <<
"#########\n";
1298 stream <<
" #########" << padded_secondary <<
"#########\n";
1299 stream <<
" ####################################################\n";
1300 stream << std::endl;
1305 template <
typename Description,
int dim,
typename Number>
1309 unsigned int timer_cycle,
1310 bool write_to_logfile,
1313 static const std::string vectorization_name = [] {
1314 constexpr auto width = VectorizedArray<Number>::size();
1320 result = std::to_string(width * 8 *
sizeof(Number)) +
" bit packed ";
1322 if constexpr (std::is_same_v<Number, double>)
1323 return result +
"double";
1324 else if constexpr (std::is_same_v<Number, float>)
1325 return result +
"float";
1330 std::ostringstream output;
1332 std::ostringstream primary;
1334 primary <<
"FINAL (cycle " << Utilities::int_to_string(cycle, 6) <<
")";
1336 primary <<
"Cycle " << Utilities::int_to_string(cycle, 6)
1337 <<
" (" << std::fixed << std::setprecision(1)
1338 << t / t_final_ * 100 <<
"%)";
1341 std::ostringstream secondary;
1342 secondary <<
"at time t = " << std::setprecision(8) << std::fixed << t;
1344 print_head(primary.str(), secondary.str(), output);
1346 output <<
"Information: (HYP) " << hyperbolic_system_.get().problem_name;
1347 if constexpr (!ParabolicSystem::is_identity) {
1348 output <<
"\n (PAR) " << parabolic_system_.get().problem_name;
1350 output <<
"\n [" << base_name_ <<
"] ";
1351 if (mpi_ensemble_.n_ensembles() > 1) {
1352 output << mpi_ensemble_.n_ensembles() <<
" ensembles ";
1355 << n_global_dofs_ <<
" Qdofs on "
1356 << mpi_ensemble_.n_world_ranks() <<
" ranks / "
1358 << MultithreadInfo::n_threads() <<
" threads <"
1360 <<
"[openmp disabled] <"
1362 << vectorization_name
1363 <<
">\n Last output cycle "
1365 <<
" at t = " << timer_granularity_ * (timer_cycle - 1)
1366 <<
" (terminal update interval " << terminal_update_interval_
1369 print_memory_statistics(output);
1370 print_timers(output);
1371 print_throughput(cycle, t, output, final_time);
1373 if (mpi_ensemble_.world_rank() == 0) {
1375 std::cout <<
"\033[2J\033[H";
1377 std::cout << output.str() << std::flush;
1379 if (write_to_logfile) {
1380 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)
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_cycle_statistics(unsigned int cycle, Number t, unsigned int output_cycle, bool write_to_logfile=false, 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)