15#include <deal.II/lac/linear_operator.h>
16#include <deal.II/lac/precondition.h>
17#include <deal.II/lac/solver_cg.h>
18#include <deal.II/matrix_free/fe_evaluation.h>
19#include <deal.II/multigrid/mg_coarse.h>
20#include <deal.II/multigrid/mg_matrix.h>
21#include <deal.II/multigrid/mg_transfer.templates.h>
22#include <deal.II/multigrid/mg_transfer_matrix_free.h>
23#include <deal.II/multigrid/multigrid.h>
29 namespace NavierStokes
31 using namespace dealii;
33 template <
typename Description,
int dim,
typename Number>
36 std::map<std::string, dealii::Timer> &computing_timer,
41 const std::string &subsection )
42 : ParameterAcceptor(subsection)
43 , mpi_ensemble_(mpi_ensemble)
44 , computing_timer_(computing_timer)
45 , hyperbolic_system_(&hyperbolic_system)
46 , parabolic_system_(¶bolic_system)
47 , offline_data_(&offline_data)
48 , initial_values_(&initial_values)
52 , n_iterations_velocity_(0.)
53 , n_iterations_internal_energy_(0.)
55 use_gmg_velocity_ =
false;
56 add_parameter(
"multigrid velocity",
58 "Use geometric multigrid for velocity component");
60 gmg_max_iter_vel_ = 12;
61 add_parameter(
"multigrid velocity - max iter",
63 "Maximal number of CG iterations with GMG smoother");
65 gmg_smoother_range_vel_ = 8.;
66 add_parameter(
"multigrid velocity - chebyshev range",
67 gmg_smoother_range_vel_,
68 "Chebyshev smoother: eigenvalue range parameter");
70 gmg_smoother_max_eig_vel_ = 2.0;
71 add_parameter(
"multigrid velocity - chebyshev max eig",
72 gmg_smoother_max_eig_vel_,
73 "Chebyshev smoother: maximal eigenvalue");
75 use_gmg_internal_energy_ =
false;
76 add_parameter(
"multigrid energy",
77 use_gmg_internal_energy_,
78 "Use geometric multigrid for internal energy component");
80 gmg_max_iter_en_ = 15;
81 add_parameter(
"multigrid energy - max iter",
83 "Maximal number of CG iterations with GMG smoother");
85 gmg_smoother_range_en_ = 15.;
86 add_parameter(
"multigrid energy - chebyshev range",
87 gmg_smoother_range_en_,
88 "Chebyshev smoother: eigenvalue range parameter");
90 gmg_smoother_max_eig_en_ = 2.0;
91 add_parameter(
"multigrid energy - chebyshev max eig",
92 gmg_smoother_max_eig_en_,
93 "Chebyshev smoother: maximal eigenvalue");
95 gmg_smoother_degree_ = 3;
96 add_parameter(
"multigrid - chebyshev degree",
98 "Chebyshev smoother: degree");
100 gmg_smoother_n_cg_iter_ = 10;
102 "multigrid - chebyshev cg iter",
103 gmg_smoother_n_cg_iter_,
104 "Chebyshev smoother: number of CG iterations to approximate "
109 "multigrid - min level",
111 "Minimal mesh level to be visited in the geometric multigrid "
112 "cycle where the coarse grid solver (Chebyshev) is called");
114 tolerance_ = Number(1.0e-12);
115 add_parameter(
"tolerance", tolerance_,
"Tolerance for linear solvers");
117 tolerance_linfty_norm_ =
false;
118 add_parameter(
"tolerance linfty norm",
119 tolerance_linfty_norm_,
120 "Use the l_infty norm instead of the l_2 norm for the "
121 "stopping criterion");
125 template <
typename Description,
int dim,
typename Number>
129 std::cout <<
"ParabolicSolver<dim, Number>::prepare()" << std::endl;
132 const auto &discretization = offline_data_->discretization();
134 dealii::ExcMessage(
"The Navier-Stokes module currently only "
135 "supports cG Q1 finite elements."));
137 AssertThrow(!offline_data_->dof_handler().has_hp_capabilities(),
139 "The Navier-Stokes module currently does not support "
140 "DofHandlers set up with hp capabilities."));
144 typename MatrixFree<dim, Number>::AdditionalData additional_data;
145 additional_data.tasks_parallel_scheme =
146 MatrixFree<dim, Number>::AdditionalData::none;
148 matrix_free_.reinit(discretization.mapping(),
149 offline_data_->dof_handler(),
150 offline_data_->affine_constraints(),
151 discretization.quadrature_1d(),
154 const auto &scalar_partitioner =
155 matrix_free_.get_dof_info(0).vector_partitioner;
157 velocity_.reinit(dim);
158 velocity_rhs_.reinit(dim);
159 for (
unsigned int i = 0; i < dim; ++i) {
160 velocity_.block(i).reinit(scalar_partitioner);
161 velocity_rhs_.block(i).reinit(scalar_partitioner);
164 internal_energy_.reinit(scalar_partitioner);
165 internal_energy_rhs_.reinit(scalar_partitioner);
167 density_.reinit(scalar_partitioner);
171 if (!use_gmg_velocity_ && !use_gmg_internal_energy_)
174 const unsigned int n_levels =
175 offline_data_->dof_handler().get_triangulation().n_global_levels();
176 const unsigned int min_level = std::min(gmg_min_level_, n_levels - 1);
177 MGLevelObject<IndexSet> relevant_sets(0, n_levels - 1);
178 for (
unsigned int level = 0; level < n_levels; ++level)
179#
if DEAL_II_VERSION_GTE(9, 6, 0)
180 relevant_sets[level] =
181 dealii::DoFTools::extract_locally_relevant_level_dofs(
182 offline_data_->dof_handler(), level);
184 dealii::DoFTools::extract_locally_relevant_level_dofs(
185 offline_data_->dof_handler(), level, relevant_sets[level]);
187 mg_constrained_dofs_.initialize(offline_data_->dof_handler(),
189 std::set<types::boundary_id> boundary_ids;
192 mg_constrained_dofs_.make_zero_boundary_constraints(
193 offline_data_->dof_handler(), boundary_ids);
195 typename MatrixFree<dim, float>::AdditionalData additional_data_level;
196 additional_data_level.tasks_parallel_scheme =
197 MatrixFree<dim, float>::AdditionalData::none;
199 level_matrix_free_.resize(min_level, n_levels - 1);
200 level_density_.resize(min_level, n_levels - 1);
201 for (
unsigned int level = min_level; level < n_levels; ++level) {
202 additional_data_level.mg_level = level;
203#if DEAL_II_VERSION_GTE(9, 6, 0)
204 AffineConstraints<double> constraints(relevant_sets[level],
205 relevant_sets[level]);
207 AffineConstraints<double> constraints(relevant_sets[level]);
212 level_matrix_free_[level].reinit(discretization.mapping(),
213 offline_data_->dof_handler(),
215 discretization.quadrature_1d(),
216 additional_data_level);
217 level_matrix_free_[level].initialize_dof_vector(level_density_[level]);
220 mg_transfer_velocity_.build(offline_data_->dof_handler(),
221 mg_constrained_dofs_,
223 mg_transfer_energy_.build(offline_data_->dof_handler(),
227 template <
typename Description,
int dim,
typename Number>
238 template <
typename Description,
int dim,
typename Number>
245 const bool reinitialize_gmg)
const
249 step(old_state_vector,
253 id_violation_strategy,
259 template <
typename Description,
int dim,
typename Number>
266 const bool reinitialize_gmg)
const
269 step(old_state_vector,
273 id_violation_strategy,
287 step(old_state_vector,
291 id_violation_strategy,
298 template <
typename Description,
int dim,
typename Number>
305 const bool reinitialize_gmg,
306 const bool crank_nicolson_extrapolation)
const
309 std::cout <<
"ParabolicSolver<dim, Number>::step()" << std::endl;
312 constexpr ScalarNumber eps = std::numeric_limits<ScalarNumber>::epsilon();
314 const auto &old_U = std::get<0>(old_state_vector);
315 auto &new_U = std::get<0>(new_state_vector);
319 using VA = VectorizedArray<Number>;
321 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
322 const auto &affine_constraints = offline_data_->affine_constraints();
326 constexpr auto simd_length = VA::size();
327 const unsigned int n_owned = offline_data_->n_locally_owned();
328 const unsigned int n_regular = n_owned / simd_length * simd_length;
330 const auto &sparsity_simd = offline_data_->sparsity_pattern_simd();
335 std::cout <<
" perform time-step with tau = " << tau << std::endl;
336 if (crank_nicolson_extrapolation)
337 std::cout <<
" and extrapolate to t + 2 * tau" << std::endl;
353 std::atomic<bool> restart_needed =
false;
364 std::atomic<bool> correction_needed =
false;
374 Scope scope(computing_timer_,
"time step [P] 1 - update velocities");
378 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
379 using T =
decltype(sentinel);
380 unsigned int stride_size = get_stride_size<T>;
382 const auto view = hyperbolic_system_->template view<dim, T>();
385 for (
unsigned int i = left; i < right; i += stride_size) {
386 const auto U_i = old_U.template get_tensor<T>(i);
387 const auto rho_i = view.density(U_i);
388 const auto M_i = view.momentum(U_i);
389 const auto rho_e_i = view.internal_energy(U_i);
390 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
392 write_entry<T>(density_, rho_i, i);
394 for (
unsigned int d = 0; d < dim; ++d) {
395 write_entry<T>(velocity_.block(d), M_i[d] / rho_i, i);
396 write_entry<T>(velocity_rhs_.block(d), m_i * (M_i[d]), i);
398 write_entry<T>(internal_energy_, rho_e_i / rho_i, i);
403 loop(Number(), n_regular, n_owned);
405 loop(VA(), 0, n_regular);
416 const auto &boundary_map = offline_data_->boundary_map();
418 for (
auto entry : boundary_map) {
420 const auto i = std::get<0>(entry);
424 const auto normal = std::get<1>(entry);
425 const auto id = std::get<4>(entry);
426 const auto position = std::get<5>(entry);
430 Tensor<1, dim, Number> V_i;
431 Tensor<1, dim, Number> RHS_i;
432 for (
unsigned int d = 0; d < dim; ++d) {
433 V_i[d] = velocity_.block(d).local_element(i);
434 RHS_i[d] = velocity_rhs_.block(d).local_element(i);
436 V_i -= 1. * (V_i * normal) * normal;
437 RHS_i -= 1. * (RHS_i * normal) * normal;
438 for (
unsigned int d = 0; d < dim; ++d) {
439 velocity_.block(d).local_element(i) = V_i[d];
440 velocity_rhs_.block(d).local_element(i) = RHS_i[d];
446 for (
unsigned int d = 0; d < dim; ++d) {
447 velocity_.block(d).local_element(i) = Number(0.);
448 velocity_rhs_.block(d).local_element(i) = Number(0.);
454 const auto U_i = initial_values_->initial_state(position, t + tau);
455 const auto view = hyperbolic_system_->template view<dim, Number>();
456 const auto rho_i = view.density(U_i);
457 const auto V_i = view.momentum(U_i) / rho_i;
458 const auto e_i = view.internal_energy(U_i) / rho_i;
460 for (
unsigned int d = 0; d < dim; ++d) {
461 velocity_.block(d).local_element(i) = V_i[d];
462 velocity_rhs_.block(d).local_element(i) = V_i[d];
465 internal_energy_.local_element(i) = e_i;
476 affine_constraints.set_zero(density_);
477 affine_constraints.set_zero(internal_energy_);
478 for (
unsigned int d = 0; d < dim; ++d) {
479 affine_constraints.set_zero(velocity_.block(d));
480 affine_constraints.set_zero(velocity_rhs_.block(d));
486 lumped_mass_matrix, density_, affine_constraints);
493 if (use_gmg_velocity_ && reinitialize_gmg) {
494 MGLevelObject<
typename PreconditionChebyshev<
495 VelocityMatrix<dim, float, Number>,
496 LinearAlgebra::distributed::BlockVector<float>,
497 DiagonalMatrix<dim, float>>::AdditionalData>
498 smoother_data(level_matrix_free_.min_level(),
499 level_matrix_free_.max_level());
501 level_velocity_matrices_.resize(level_matrix_free_.min_level(),
502 level_matrix_free_.max_level());
503 mg_transfer_velocity_.interpolate_to_mg(
504 offline_data_->dof_handler(), level_density_, density_);
506 for (
unsigned int level = level_matrix_free_.min_level();
507 level <= level_matrix_free_.max_level();
509 level_velocity_matrices_[level].initialize(
512 level_matrix_free_[level],
513 level_density_[level],
516 level_velocity_matrices_[level].compute_diagonal(
517 smoother_data[level].preconditioner);
518 if (level == level_matrix_free_.min_level()) {
519 smoother_data[level].degree = numbers::invalid_unsigned_int;
520 smoother_data[level].eig_cg_n_iterations = 500;
521 smoother_data[level].smoothing_range = 1e-3;
523 smoother_data[level].degree = gmg_smoother_degree_;
524 smoother_data[level].eig_cg_n_iterations =
525 gmg_smoother_n_cg_iter_;
526 smoother_data[level].smoothing_range = gmg_smoother_range_vel_;
527 if (gmg_smoother_n_cg_iter_ == 0)
528 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_vel_;
531 mg_smoother_velocity_.initialize(level_velocity_matrices_,
541 Scope scope(computing_timer_,
542 "time step [X] _ - synchronization barriers");
548 *std::min_element(internal_energy_.begin(), internal_energy_.end());
550 e_min_old = Utilities::MPI::min(e_min_old,
551 mpi_ensemble_.ensemble_communicator());
554 constexpr Number eps = std::numeric_limits<Number>::epsilon();
555 e_min_old *= (1. - 1000. * eps);
562 Scope scope(computing_timer_,
"time step [P] 1 - update velocities");
566 VelocityMatrix<dim, Number, Number> velocity_operator;
567 velocity_operator.initialize(
568 *parabolic_system_, *offline_data_, matrix_free_, density_, tau);
570 const auto tolerance_velocity =
571 (tolerance_linfty_norm_ ? velocity_rhs_.linfty_norm()
572 : velocity_rhs_.l2_norm()) *
581 if (!use_gmg_velocity_)
582 throw SolverControl::NoConvergence(0, 0.);
584 using bvt_float = LinearAlgebra::distributed::BlockVector<float>;
586 MGCoarseGridApplySmoother<bvt_float> mg_coarse;
587 mg_coarse.initialize(mg_smoother_velocity_);
589 mg::Matrix<bvt_float> mg_matrix(level_velocity_matrices_);
591 Multigrid<bvt_float> mg(mg_matrix,
593 mg_transfer_velocity_,
594 mg_smoother_velocity_,
595 mg_smoother_velocity_,
596 level_velocity_matrices_.min_level(),
597 level_velocity_matrices_.max_level());
599 const auto &dof_handler = offline_data_->dof_handler();
600 PreconditionMG<dim, bvt_float, MGTransferVelocity<dim, float>>
601 preconditioner(dof_handler, mg, mg_transfer_velocity_);
603 SolverControl solver_control(gmg_max_iter_vel_, tolerance_velocity);
604 SolverCG<BlockVector> solver(solver_control);
606 velocity_operator, velocity_, velocity_rhs_, preconditioner);
609 n_iterations_velocity_ =
610 0.9 * n_iterations_velocity_ + 0.1 * solver_control.last_step();
612 }
catch (SolverControl::NoConvergence &) {
614 SolverControl solver_control(1000, tolerance_velocity);
615 SolverCG<BlockVector> solver(solver_control);
617 velocity_operator, velocity_, velocity_rhs_, diagonal_matrix);
620 n_iterations_velocity_ *= 0.9;
621 n_iterations_velocity_ +=
622 0.1 * (use_gmg_velocity_ ? gmg_max_iter_vel_ : 0) +
623 0.1 * solver_control.last_step();
633 Scope scope(computing_timer_,
634 "time step [P] 2 - update internal energy");
639 matrix_free_.template cell_loop<ScalarVector, BlockVector>(
640 [
this](
const auto &data,
643 const auto cell_range) {
644 FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
646 FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(data);
648 const auto mu = parabolic_system_->mu();
649 const auto lambda = parabolic_system_->lambda();
651 for (
unsigned int cell = cell_range.first;
652 cell < cell_range.second;
654 velocity.reinit(cell);
656 velocity.gather_evaluate(src, EvaluationFlags::gradients);
658 for (
unsigned int q = 0; q < velocity.n_q_points; ++q) {
659 if constexpr (dim == 1) {
661 const auto gradient = velocity.get_gradient(q);
662 auto S = (4. / 3. * mu + lambda) * gradient;
663 energy.submit_value(gradient * S, q);
667 const auto symmetric_gradient =
668 velocity.get_symmetric_gradient(q);
669 const auto divergence = trace(symmetric_gradient);
670 auto S = 2. * mu * symmetric_gradient;
671 for (
unsigned int d = 0; d < dim; ++d)
672 S[d][d] += (lambda - 2. / 3. * mu) * divergence;
673 energy.submit_value(symmetric_gradient * S, q);
676 energy.integrate_scatter(EvaluationFlags::values, dst);
679 internal_energy_rhs_,
683 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
687 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
688 using T =
decltype(sentinel);
689 unsigned int stride_size = get_stride_size<T>;
691 const auto view = hyperbolic_system_->template view<dim, T>();
694 for (
unsigned int i = left; i < right; i += stride_size) {
695 const auto rhs_i = get_entry<T>(internal_energy_rhs_, i);
696 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
697 const auto rho_i = get_entry<T>(density_, i);
698 const auto e_i = get_entry<T>(internal_energy_, i);
700 const auto U_i = old_U.template get_tensor<T>(i);
701 const auto V_i = view.momentum(U_i) / rho_i;
703 dealii::Tensor<1, dim, T> V_i_new;
704 for (
unsigned int d = 0; d < dim; ++d) {
705 V_i_new[d] = get_entry<T>(velocity_.block(d), i);
712 const auto correction =
713 crank_nicolson_extrapolation
715 : Number(0.5) * (V_i - V_i_new).norm_square();
718 const auto result = m_i * rho_i * (e_i + correction) + tau * rhs_i;
719 write_entry<T>(internal_energy_rhs_, result, i);
724 loop(Number(), n_regular, n_owned);
726 loop(VA(), 0, n_regular);
737 const auto &boundary_map = offline_data_->boundary_map();
739 for (
auto entry : boundary_map) {
741 const auto i = std::get<0>(entry);
745 const auto id = std::get<4>(entry);
746 const auto position = std::get<5>(entry);
750 const auto U_i = initial_values_->initial_state(position, t + tau);
751 const auto view = hyperbolic_system_->template view<dim, Number>();
752 const auto rho_i = view.density(U_i);
753 const auto e_i = view.internal_energy(U_i) / rho_i;
754 internal_energy_rhs_.local_element(i) = e_i;
764 affine_constraints.set_zero(internal_energy_);
765 affine_constraints.set_zero(internal_energy_rhs_);
772 if (use_gmg_internal_energy_ && reinitialize_gmg) {
773 MGLevelObject<
typename PreconditionChebyshev<
774 EnergyMatrix<dim, float, Number>,
775 LinearAlgebra::distributed::Vector<float>>::AdditionalData>
776 smoother_data(level_matrix_free_.min_level(),
777 level_matrix_free_.max_level());
779 level_energy_matrices_.resize(level_matrix_free_.min_level(),
780 level_matrix_free_.max_level());
782 for (
unsigned int level = level_matrix_free_.min_level();
783 level <= level_matrix_free_.max_level();
785 level_energy_matrices_[level].initialize(
787 level_matrix_free_[level],
788 level_density_[level],
789 tau * parabolic_system_->cv_inverse_kappa(),
791 level_energy_matrices_[level].compute_diagonal(
792 smoother_data[level].preconditioner);
793 if (level == level_matrix_free_.min_level()) {
794 smoother_data[level].degree = numbers::invalid_unsigned_int;
795 smoother_data[level].eig_cg_n_iterations = 500;
796 smoother_data[level].smoothing_range = 1e-3;
798 smoother_data[level].degree = gmg_smoother_degree_;
799 smoother_data[level].eig_cg_n_iterations =
800 gmg_smoother_n_cg_iter_;
801 smoother_data[level].smoothing_range = gmg_smoother_range_en_;
802 if (gmg_smoother_n_cg_iter_ == 0)
803 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_en_;
806 mg_smoother_energy_.initialize(level_energy_matrices_, smoother_data);
816 Scope scope(computing_timer_,
817 "time step [P] 2 - update internal energy");
821 EnergyMatrix<dim, Number, Number> energy_operator;
822 const auto &kappa = parabolic_system_->cv_inverse_kappa();
823 energy_operator.initialize(
824 *offline_data_, matrix_free_, density_, tau * kappa);
826 const auto tolerance_internal_energy =
827 (tolerance_linfty_norm_ ? internal_energy_rhs_.linfty_norm()
828 : internal_energy_rhs_.l2_norm()) *
832 if (!use_gmg_internal_energy_)
833 throw SolverControl::NoConvergence(0, 0.);
835 using vt_float = LinearAlgebra::distributed::Vector<float>;
836 MGCoarseGridApplySmoother<vt_float> mg_coarse;
837 mg_coarse.initialize(mg_smoother_energy_);
838 mg::Matrix<vt_float> mg_matrix(level_energy_matrices_);
840 Multigrid<vt_float> mg(mg_matrix,
845 level_energy_matrices_.min_level(),
846 level_energy_matrices_.max_level());
848 const auto &dof_handler = offline_data_->dof_handler();
849 PreconditionMG<dim, vt_float, MGTransferEnergy<dim, float>>
850 preconditioner(dof_handler, mg, mg_transfer_energy_);
852 SolverControl solver_control(gmg_max_iter_en_,
853 tolerance_internal_energy);
854 SolverCG<ScalarVector> solver(solver_control);
855 solver.solve(energy_operator,
857 internal_energy_rhs_,
861 n_iterations_internal_energy_ = 0.9 * n_iterations_internal_energy_ +
862 0.1 * solver_control.last_step();
864 }
catch (SolverControl::NoConvergence &) {
866 SolverControl solver_control(1000, tolerance_internal_energy);
867 SolverCG<ScalarVector> solver(solver_control);
868 solver.solve(energy_operator,
870 internal_energy_rhs_,
874 n_iterations_internal_energy_ *= 0.9;
875 n_iterations_internal_energy_ +=
876 0.1 * (use_gmg_internal_energy_ ? gmg_max_iter_en_ : 0) +
877 0.1 * solver_control.last_step();
890 Scope scope(computing_timer_,
"time step [P] 3 - write back vectors");
895 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
896 using T =
decltype(sentinel);
897 unsigned int stride_size = get_stride_size<T>;
899 const auto view = hyperbolic_system_->template view<dim, T>();
902 for (
unsigned int i = left; i < right; i += stride_size) {
905 const unsigned int row_length = sparsity_simd.row_length(i);
909 auto U_i = old_U.template get_tensor<T>(i);
910 const auto rho_i = view.density(U_i);
912 Tensor<1, dim, T> m_i_new;
913 for (
unsigned int d = 0; d < dim; ++d) {
914 m_i_new[d] = rho_i * get_entry<T>(velocity_.block(d), i);
917 auto rho_e_i_new = rho_i * get_entry<T>(internal_energy_, i);
924 if (!(T(0.) == std::max(T(0.), rho_i * e_min_old - rho_e_i_new))) {
926 std::cout << std::fixed << std::setprecision(16);
927 const auto e_i_new = rho_e_i_new / rho_i;
928 std::cout <<
"Bounds violation: internal energy (critical)!\n"
929 <<
"\t\te_min_old: " << e_min_old <<
"\n"
930 <<
"\t\te_min_old (delta): "
932 <<
"\t\te_min_new: " << e_i_new <<
"\n"
935 restart_needed =
true;
938 if (crank_nicolson_extrapolation) {
939 m_i_new = Number(2.0) * m_i_new - view.momentum(U_i);
941 Number(2.0) * rho_e_i_new - view.internal_energy(U_i);
949 std::max(T(0.), eps * rho_i * e_min_old - rho_e_i_new))) {
951 std::cout << std::fixed << std::setprecision(16);
952 const auto e_i_new = rho_e_i_new / rho_i;
954 std::cout <<
"Bounds violation: high-order internal energy!"
955 <<
"\t\te_min_new: " << e_i_new <<
"\n"
956 <<
"\t\t-- correction required --" << std::endl;
958 correction_needed =
true;
962 const auto E_i_new = rho_e_i_new + 0.5 * m_i_new * m_i_new / rho_i;
964 for (
unsigned int d = 0; d < dim; ++d)
965 U_i[1 + d] = m_i_new[d];
966 U_i[1 + dim] = E_i_new;
968 new_U.template write_tensor<T>(U_i, i);
973 loop(Number(), n_regular, n_owned);
975 loop(VA(), 0, n_regular);
980 new_U.update_ghost_values();
986 Scope scope(computing_timer_,
987 "time step [X] _ - synchronization barriers");
998 restart_needed.store(Utilities::MPI::logical_or(
999 restart_needed.load(),
1000 mpi_ensemble_.synchronization_communicator()));
1002 correction_needed.store(Utilities::MPI::logical_or(
1003 correction_needed.load(),
1004 mpi_ensemble_.synchronization_communicator()));
1007 if (correction_needed) {
1012 throw Restart{Number(0.5) * tau};
1019 if (restart_needed) {
1020 switch (id_violation_strategy) {
1027 throw Restart{Number(0.5) * tau};
1033 template <
typename Description,
int dim,
typename Number>
1035 std::ostream &output)
const
1037 output <<
" [ " << std::setprecision(2) << std::fixed
1038 << n_iterations_velocity_
1039 << (use_gmg_velocity_ ?
" GMG vel -- " :
" CG vel -- ")
1040 << n_iterations_internal_energy_
1041 << (use_gmg_internal_energy_ ?
" GMG int ]" :
" CG int ]")
void reinit(const vector_type &lumped_mass_matrix, const vector_type &density, const dealii::AffineConstraints< Number > &affine_constraints)
void backward_euler_step(const StateVector &old_state_vector, const Number old_t, StateVector &new_state_vector, Number tau, const IDViolationStrategy id_violation_strategy, const bool reinitialize_gmg) const
typename Description::ParabolicSystem ParabolicSystem
typename View::StateVector StateVector
void print_solver_statistics(std::ostream &output) const
void prepare_state_vector(StateVector &state_vector, Number t) const
typename Description::HyperbolicSystem HyperbolicSystem
ParabolicSolver(const MPIEnsemble &mpi_ensemble, std::map< std::string, dealii::Timer > &computing_timer, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem ¶bolic_system, const OfflineData< dim, Number > &offline_data, const InitialValues< Description, dim, Number > &initial_values, const std::string &subsection="ParabolicSolver")
void crank_nicolson_step(const StateVector &old_state_vector, const Number old_t, StateVector &new_state_vector, Number tau, const IDViolationStrategy id_violation_strategy, const bool reinitialize_gmg) const
void step(Triangulation< dim, dim > &, const double, const double, const double, const double)
#define RYUJIN_PARALLEL_REGION_BEGIN
#define RYUJIN_PARALLEL_REGION_END
DEAL_II_ALWAYS_INLINE Number negative_part(const Number number)
#define LIKWID_MARKER_START(opt)
#define CALLGRIND_START_INSTRUMENTATION
#define LIKWID_MARKER_STOP(opt)
#define CALLGRIND_STOP_INSTRUMENTATION
std::tuple< MultiComponentVector< Number, problem_dim >, MultiComponentVector< Number, prec_dim >, BlockVector< Number > > StateVector