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 NavierStokes module currently only "
135 "supports cG Q1 finite elements."));
139 typename MatrixFree<dim, Number>::AdditionalData additional_data;
140 additional_data.tasks_parallel_scheme =
141 MatrixFree<dim, Number>::AdditionalData::none;
143 matrix_free_.reinit(discretization.mapping(),
144 offline_data_->dof_handler(),
145 offline_data_->affine_constraints(),
146 discretization.quadrature_1d(),
149 const auto &scalar_partitioner =
150 matrix_free_.get_dof_info(0).vector_partitioner;
152 velocity_.reinit(dim);
153 velocity_rhs_.reinit(dim);
154 for (
unsigned int i = 0; i < dim; ++i) {
155 velocity_.block(i).reinit(scalar_partitioner);
156 velocity_rhs_.block(i).reinit(scalar_partitioner);
159 internal_energy_.reinit(scalar_partitioner);
160 internal_energy_rhs_.reinit(scalar_partitioner);
162 density_.reinit(scalar_partitioner);
166 if (!use_gmg_velocity_ && !use_gmg_internal_energy_)
169 const unsigned int n_levels =
170 offline_data_->dof_handler().get_triangulation().n_global_levels();
171 const unsigned int min_level = std::min(gmg_min_level_, n_levels - 1);
172 MGLevelObject<IndexSet> relevant_sets(0, n_levels - 1);
173 for (
unsigned int level = 0; level < n_levels; ++level)
174 dealii::DoFTools::extract_locally_relevant_level_dofs(
175 offline_data_->dof_handler(), level, relevant_sets[level]);
176 mg_constrained_dofs_.initialize(offline_data_->dof_handler(),
178 std::set<types::boundary_id> boundary_ids;
181 mg_constrained_dofs_.make_zero_boundary_constraints(
182 offline_data_->dof_handler(), boundary_ids);
184 typename MatrixFree<dim, float>::AdditionalData additional_data_level;
185 additional_data_level.tasks_parallel_scheme =
186 MatrixFree<dim, float>::AdditionalData::none;
188 level_matrix_free_.resize(min_level, n_levels - 1);
189 level_density_.resize(min_level, n_levels - 1);
190 for (
unsigned int level = min_level; level < n_levels; ++level) {
191 additional_data_level.mg_level = level;
192 AffineConstraints<double> constraints(relevant_sets[level]);
196 level_matrix_free_[level].reinit(discretization.mapping(),
197 offline_data_->dof_handler(),
199 discretization.quadrature_1d(),
200 additional_data_level);
201 level_matrix_free_[level].initialize_dof_vector(level_density_[level]);
204 mg_transfer_velocity_.build(offline_data_->dof_handler(),
205 mg_constrained_dofs_,
207 mg_transfer_energy_.build(offline_data_->dof_handler(),
211 template <
typename Description,
int dim,
typename Number>
222 template <
typename Description,
int dim,
typename Number>
229 const bool reinitialize_gmg)
const
233 step(old_state_vector,
237 id_violation_strategy,
243 template <
typename Description,
int dim,
typename Number>
250 const bool reinitialize_gmg)
const
253 step(old_state_vector,
257 id_violation_strategy,
271 step(old_state_vector,
275 id_violation_strategy,
282 template <
typename Description,
int dim,
typename Number>
289 const bool reinitialize_gmg,
290 const bool crank_nicolson_extrapolation)
const
293 std::cout <<
"ParabolicSolver<dim, Number>::step()" << std::endl;
296 constexpr ScalarNumber eps = std::numeric_limits<ScalarNumber>::epsilon();
298 const auto &old_U = std::get<0>(old_state_vector);
299 auto &new_U = std::get<0>(new_state_vector);
303 using VA = VectorizedArray<Number>;
305 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
306 const auto &affine_constraints = offline_data_->affine_constraints();
310 constexpr auto simd_length = VA::size();
311 const unsigned int n_owned = offline_data_->n_locally_owned();
312 const unsigned int n_regular = n_owned / simd_length * simd_length;
314 const auto &sparsity_simd = offline_data_->sparsity_pattern_simd();
319 std::cout <<
" perform time-step with tau = " << tau << std::endl;
320 if (crank_nicolson_extrapolation)
321 std::cout <<
" and extrapolate to t + 2 * tau" << std::endl;
337 std::atomic<bool> restart_needed =
false;
348 std::atomic<bool> correction_needed =
false;
358 Scope scope(computing_timer_,
"time step [P] 1 - update velocities");
362 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
363 using T =
decltype(sentinel);
364 unsigned int stride_size = get_stride_size<T>;
366 const auto view = hyperbolic_system_->template view<dim, T>();
369 for (
unsigned int i = left; i < right; i += stride_size) {
370 const auto U_i = old_U.template get_tensor<T>(i);
371 const auto rho_i = view.density(U_i);
372 const auto M_i = view.momentum(U_i);
373 const auto rho_e_i = view.internal_energy(U_i);
374 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
376 write_entry<T>(density_, rho_i, i);
378 for (
unsigned int d = 0; d < dim; ++d) {
379 write_entry<T>(velocity_.block(d), M_i[d] / rho_i, i);
380 write_entry<T>(velocity_rhs_.block(d), m_i * (M_i[d]), i);
382 write_entry<T>(internal_energy_, rho_e_i / rho_i, i);
387 loop(Number(), n_regular, n_owned);
389 loop(VA(), 0, n_regular);
400 const auto &boundary_map = offline_data_->boundary_map();
402 for (
auto entry : boundary_map) {
404 const auto i = std::get<0>(entry);
408 const auto normal = std::get<1>(entry);
409 const auto id = std::get<4>(entry);
410 const auto position = std::get<5>(entry);
414 Tensor<1, dim, Number> V_i;
415 Tensor<1, dim, Number> RHS_i;
416 for (
unsigned int d = 0; d < dim; ++d) {
417 V_i[d] = velocity_.block(d).local_element(i);
418 RHS_i[d] = velocity_rhs_.block(d).local_element(i);
420 V_i -= 1. * (V_i * normal) * normal;
421 RHS_i -= 1. * (RHS_i * normal) * normal;
422 for (
unsigned int d = 0; d < dim; ++d) {
423 velocity_.block(d).local_element(i) = V_i[d];
424 velocity_rhs_.block(d).local_element(i) = RHS_i[d];
430 for (
unsigned int d = 0; d < dim; ++d) {
431 velocity_.block(d).local_element(i) = Number(0.);
432 velocity_rhs_.block(d).local_element(i) = Number(0.);
438 const auto U_i = initial_values_->initial_state(position, t + tau);
439 const auto view = hyperbolic_system_->template view<dim, Number>();
440 const auto rho_i = view.density(U_i);
441 const auto V_i = view.momentum(U_i) / rho_i;
442 const auto e_i = view.internal_energy(U_i) / rho_i;
444 for (
unsigned int d = 0; d < dim; ++d) {
445 velocity_.block(d).local_element(i) = V_i[d];
446 velocity_rhs_.block(d).local_element(i) = V_i[d];
449 internal_energy_.local_element(i) = e_i;
460 affine_constraints.set_zero(density_);
461 affine_constraints.set_zero(internal_energy_);
462 for (
unsigned int d = 0; d < dim; ++d) {
463 affine_constraints.set_zero(velocity_.block(d));
464 affine_constraints.set_zero(velocity_rhs_.block(d));
470 lumped_mass_matrix, density_, affine_constraints);
477 if (use_gmg_velocity_ && reinitialize_gmg) {
478 MGLevelObject<
typename PreconditionChebyshev<
479 VelocityMatrix<dim, float, Number>,
480 LinearAlgebra::distributed::BlockVector<float>,
481 DiagonalMatrix<dim, float>>::AdditionalData>
482 smoother_data(level_matrix_free_.min_level(),
483 level_matrix_free_.max_level());
485 level_velocity_matrices_.resize(level_matrix_free_.min_level(),
486 level_matrix_free_.max_level());
487 mg_transfer_velocity_.interpolate_to_mg(
488 offline_data_->dof_handler(), level_density_, density_);
490 for (
unsigned int level = level_matrix_free_.min_level();
491 level <= level_matrix_free_.max_level();
493 level_velocity_matrices_[level].initialize(
496 level_matrix_free_[level],
497 level_density_[level],
500 level_velocity_matrices_[level].compute_diagonal(
501 smoother_data[level].preconditioner);
502 if (level == level_matrix_free_.min_level()) {
503 smoother_data[level].degree = numbers::invalid_unsigned_int;
504 smoother_data[level].eig_cg_n_iterations = 500;
505 smoother_data[level].smoothing_range = 1e-3;
507 smoother_data[level].degree = gmg_smoother_degree_;
508 smoother_data[level].eig_cg_n_iterations =
509 gmg_smoother_n_cg_iter_;
510 smoother_data[level].smoothing_range = gmg_smoother_range_vel_;
511 if (gmg_smoother_n_cg_iter_ == 0)
512 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_vel_;
515 mg_smoother_velocity_.initialize(level_velocity_matrices_,
525 Scope scope(computing_timer_,
526 "time step [P] _ - synchronization barriers");
532 *std::min_element(internal_energy_.begin(), internal_energy_.end());
534 e_min_old = Utilities::MPI::min(e_min_old,
535 mpi_ensemble_.ensemble_communicator());
538 constexpr Number eps = std::numeric_limits<Number>::epsilon();
539 e_min_old *= (1. - 1000. * eps);
546 Scope scope(computing_timer_,
"time step [P] 1 - update velocities");
550 VelocityMatrix<dim, Number, Number> velocity_operator;
551 velocity_operator.initialize(
552 *parabolic_system_, *offline_data_, matrix_free_, density_, tau);
554 const auto tolerance_velocity =
555 (tolerance_linfty_norm_ ? velocity_rhs_.linfty_norm()
556 : velocity_rhs_.l2_norm()) *
565 if (!use_gmg_velocity_)
566 throw SolverControl::NoConvergence(0, 0.);
568 using bvt_float = LinearAlgebra::distributed::BlockVector<float>;
570 MGCoarseGridApplySmoother<bvt_float> mg_coarse;
571 mg_coarse.initialize(mg_smoother_velocity_);
573 mg::Matrix<bvt_float> mg_matrix(level_velocity_matrices_);
575 Multigrid<bvt_float> mg(mg_matrix,
577 mg_transfer_velocity_,
578 mg_smoother_velocity_,
579 mg_smoother_velocity_,
580 level_velocity_matrices_.min_level(),
581 level_velocity_matrices_.max_level());
583 const auto &dof_handler = offline_data_->dof_handler();
584 PreconditionMG<dim, bvt_float, MGTransferVelocity<dim, float>>
585 preconditioner(dof_handler, mg, mg_transfer_velocity_);
587 SolverControl solver_control(gmg_max_iter_vel_, tolerance_velocity);
588 SolverCG<BlockVector> solver(solver_control);
590 velocity_operator, velocity_, velocity_rhs_, preconditioner);
593 n_iterations_velocity_ =
594 0.9 * n_iterations_velocity_ + 0.1 * solver_control.last_step();
596 }
catch (SolverControl::NoConvergence &) {
598 SolverControl solver_control(1000, tolerance_velocity);
599 SolverCG<BlockVector> solver(solver_control);
601 velocity_operator, velocity_, velocity_rhs_, diagonal_matrix);
604 n_iterations_velocity_ *= 0.9;
605 n_iterations_velocity_ +=
606 0.1 * (use_gmg_velocity_ ? gmg_max_iter_vel_ : 0) +
607 0.1 * solver_control.last_step();
617 Scope scope(computing_timer_,
618 "time step [P] 2 - update internal energy");
623 matrix_free_.template cell_loop<ScalarVector, BlockVector>(
624 [
this](
const auto &data,
627 const auto cell_range) {
628 FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
630 FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(data);
632 const auto mu = parabolic_system_->mu();
633 const auto lambda = parabolic_system_->lambda();
635 for (
unsigned int cell = cell_range.first;
636 cell < cell_range.second;
638 velocity.reinit(cell);
640 velocity.gather_evaluate(src, EvaluationFlags::gradients);
642 for (
unsigned int q = 0; q < velocity.n_q_points; ++q) {
643 if constexpr (dim == 1) {
645 const auto gradient = velocity.get_gradient(q);
646 auto S = (4. / 3. * mu + lambda) * gradient;
647 energy.submit_value(gradient * S, q);
651 const auto symmetric_gradient =
652 velocity.get_symmetric_gradient(q);
653 const auto divergence = trace(symmetric_gradient);
654 auto S = 2. * mu * symmetric_gradient;
655 for (
unsigned int d = 0; d < dim; ++d)
656 S[d][d] += (lambda - 2. / 3. * mu) * divergence;
657 energy.submit_value(symmetric_gradient * S, q);
660 energy.integrate_scatter(EvaluationFlags::values, dst);
663 internal_energy_rhs_,
667 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
671 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
672 using T =
decltype(sentinel);
673 unsigned int stride_size = get_stride_size<T>;
675 const auto view = hyperbolic_system_->template view<dim, T>();
678 for (
unsigned int i = left; i < right; i += stride_size) {
679 const auto rhs_i = get_entry<T>(internal_energy_rhs_, i);
680 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
681 const auto rho_i = get_entry<T>(density_, i);
682 const auto e_i = get_entry<T>(internal_energy_, i);
684 const auto U_i = old_U.template get_tensor<T>(i);
685 const auto V_i = view.momentum(U_i) / rho_i;
687 dealii::Tensor<1, dim, T> V_i_new;
688 for (
unsigned int d = 0; d < dim; ++d) {
689 V_i_new[d] = get_entry<T>(velocity_.block(d), i);
696 const auto correction =
697 crank_nicolson_extrapolation
699 : Number(0.5) * (V_i - V_i_new).norm_square();
702 const auto result = m_i * rho_i * (e_i + correction) + tau * rhs_i;
703 write_entry<T>(internal_energy_rhs_, result, i);
708 loop(Number(), n_regular, n_owned);
710 loop(VA(), 0, n_regular);
721 const auto &boundary_map = offline_data_->boundary_map();
723 for (
auto entry : boundary_map) {
725 const auto i = std::get<0>(entry);
729 const auto id = std::get<4>(entry);
730 const auto position = std::get<5>(entry);
734 const auto U_i = initial_values_->initial_state(position, t + tau);
735 const auto view = hyperbolic_system_->template view<dim, Number>();
736 const auto rho_i = view.density(U_i);
737 const auto e_i = view.internal_energy(U_i) / rho_i;
738 internal_energy_rhs_.local_element(i) = e_i;
748 affine_constraints.set_zero(internal_energy_);
749 affine_constraints.set_zero(internal_energy_rhs_);
756 if (use_gmg_internal_energy_ && reinitialize_gmg) {
757 MGLevelObject<
typename PreconditionChebyshev<
758 EnergyMatrix<dim, float, Number>,
759 LinearAlgebra::distributed::Vector<float>>::AdditionalData>
760 smoother_data(level_matrix_free_.min_level(),
761 level_matrix_free_.max_level());
763 level_energy_matrices_.resize(level_matrix_free_.min_level(),
764 level_matrix_free_.max_level());
766 for (
unsigned int level = level_matrix_free_.min_level();
767 level <= level_matrix_free_.max_level();
769 level_energy_matrices_[level].initialize(
771 level_matrix_free_[level],
772 level_density_[level],
773 tau * parabolic_system_->cv_inverse_kappa(),
775 level_energy_matrices_[level].compute_diagonal(
776 smoother_data[level].preconditioner);
777 if (level == level_matrix_free_.min_level()) {
778 smoother_data[level].degree = numbers::invalid_unsigned_int;
779 smoother_data[level].eig_cg_n_iterations = 500;
780 smoother_data[level].smoothing_range = 1e-3;
782 smoother_data[level].degree = gmg_smoother_degree_;
783 smoother_data[level].eig_cg_n_iterations =
784 gmg_smoother_n_cg_iter_;
785 smoother_data[level].smoothing_range = gmg_smoother_range_en_;
786 if (gmg_smoother_n_cg_iter_ == 0)
787 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_en_;
790 mg_smoother_energy_.initialize(level_energy_matrices_, smoother_data);
800 Scope scope(computing_timer_,
801 "time step [P] 2 - update internal energy");
805 EnergyMatrix<dim, Number, Number> energy_operator;
806 const auto &kappa = parabolic_system_->cv_inverse_kappa();
807 energy_operator.initialize(
808 *offline_data_, matrix_free_, density_, tau * kappa);
810 const auto tolerance_internal_energy =
811 (tolerance_linfty_norm_ ? internal_energy_rhs_.linfty_norm()
812 : internal_energy_rhs_.l2_norm()) *
816 if (!use_gmg_internal_energy_)
817 throw SolverControl::NoConvergence(0, 0.);
819 using vt_float = LinearAlgebra::distributed::Vector<float>;
820 MGCoarseGridApplySmoother<vt_float> mg_coarse;
821 mg_coarse.initialize(mg_smoother_energy_);
822 mg::Matrix<vt_float> mg_matrix(level_energy_matrices_);
824 Multigrid<vt_float> mg(mg_matrix,
829 level_energy_matrices_.min_level(),
830 level_energy_matrices_.max_level());
832 const auto &dof_handler = offline_data_->dof_handler();
833 PreconditionMG<dim, vt_float, MGTransferEnergy<dim, float>>
834 preconditioner(dof_handler, mg, mg_transfer_energy_);
836 SolverControl solver_control(gmg_max_iter_en_,
837 tolerance_internal_energy);
838 SolverCG<ScalarVector> solver(solver_control);
839 solver.solve(energy_operator,
841 internal_energy_rhs_,
845 n_iterations_internal_energy_ = 0.9 * n_iterations_internal_energy_ +
846 0.1 * solver_control.last_step();
848 }
catch (SolverControl::NoConvergence &) {
850 SolverControl solver_control(1000, tolerance_internal_energy);
851 SolverCG<ScalarVector> solver(solver_control);
852 solver.solve(energy_operator,
854 internal_energy_rhs_,
858 n_iterations_internal_energy_ *= 0.9;
859 n_iterations_internal_energy_ +=
860 0.1 * (use_gmg_internal_energy_ ? gmg_max_iter_en_ : 0) +
861 0.1 * solver_control.last_step();
874 Scope scope(computing_timer_,
"time step [P] 3 - write back vectors");
879 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
880 using T =
decltype(sentinel);
881 unsigned int stride_size = get_stride_size<T>;
883 const auto view = hyperbolic_system_->template view<dim, T>();
886 for (
unsigned int i = left; i < right; i += stride_size) {
889 const unsigned int row_length = sparsity_simd.row_length(i);
893 auto U_i = old_U.template get_tensor<T>(i);
894 const auto rho_i = view.density(U_i);
896 Tensor<1, dim, T> m_i_new;
897 for (
unsigned int d = 0; d < dim; ++d) {
898 m_i_new[d] = rho_i * get_entry<T>(velocity_.block(d), i);
901 auto rho_e_i_new = rho_i * get_entry<T>(internal_energy_, i);
908 if (!(T(0.) == std::max(T(0.), rho_i * e_min_old - rho_e_i_new))) {
910 std::cout << std::fixed << std::setprecision(16);
911 const auto e_i_new = rho_e_i_new / rho_i;
912 std::cout <<
"Bounds violation: internal energy (critical)!\n"
913 <<
"\t\te_min_old: " << e_min_old <<
"\n"
914 <<
"\t\te_min_old (delta): "
916 <<
"\t\te_min_new: " << e_i_new <<
"\n"
919 restart_needed =
true;
922 if (crank_nicolson_extrapolation) {
923 m_i_new = Number(2.0) * m_i_new - view.momentum(U_i);
925 Number(2.0) * rho_e_i_new - view.internal_energy(U_i);
933 std::max(T(0.), eps * rho_i * e_min_old - rho_e_i_new))) {
935 std::cout << std::fixed << std::setprecision(16);
936 const auto e_i_new = rho_e_i_new / rho_i;
938 std::cout <<
"Bounds violation: high-order internal energy!"
939 <<
"\t\te_min_new: " << e_i_new <<
"\n"
940 <<
"\t\t-- correction required --" << std::endl;
942 correction_needed =
true;
946 const auto E_i_new = rho_e_i_new + 0.5 * m_i_new * m_i_new / rho_i;
948 for (
unsigned int d = 0; d < dim; ++d)
949 U_i[1 + d] = m_i_new[d];
950 U_i[1 + dim] = E_i_new;
952 new_U.template write_tensor<T>(U_i, i);
957 loop(Number(), n_regular, n_owned);
959 loop(VA(), 0, n_regular);
964 new_U.update_ghost_values();
970 Scope scope(computing_timer_,
971 "time step [H] _ - synchronization barriers");
982 restart_needed.store(Utilities::MPI::logical_or(
983 restart_needed.load(),
984 mpi_ensemble_.synchronization_communicator()));
986 correction_needed.store(Utilities::MPI::logical_or(
987 correction_needed.load(),
988 mpi_ensemble_.synchronization_communicator()));
991 if (correction_needed) {
996 throw Restart{Number(0.5) * tau};
1003 if (restart_needed) {
1004 switch (id_violation_strategy) {
1011 throw Restart{Number(0.5) * tau};
1017 template <
typename Description,
int dim,
typename Number>
1019 std::ostream &output)
const
1021 output <<
" [ " << std::setprecision(2) << std::fixed
1022 << n_iterations_velocity_
1023 << (use_gmg_velocity_ ?
" GMG vel -- " :
" CG vel -- ")
1024 << n_iterations_internal_energy_
1025 << (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