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(),
212 template <
typename Description,
int dim,
typename Number>
219 const bool reinitialize_gmg)
const
223 step(old_state_vector,
227 id_violation_strategy,
233 template <
typename Description,
int dim,
typename Number>
240 const bool reinitialize_gmg)
const
243 step(old_state_vector,
247 id_violation_strategy,
261 step(old_state_vector,
265 id_violation_strategy,
272 template <
typename Description,
int dim,
typename Number>
279 const bool reinitialize_gmg,
280 const bool crank_nicolson_extrapolation)
const
283 std::cout <<
"ParabolicSolver<dim, Number>::step()" << std::endl;
286 constexpr ScalarNumber eps = std::numeric_limits<ScalarNumber>::epsilon();
288 const auto &old_U = std::get<0>(old_state_vector);
289 auto &new_U = std::get<0>(new_state_vector);
293 using VA = VectorizedArray<Number>;
295 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
296 const auto &affine_constraints = offline_data_->affine_constraints();
300 constexpr auto simd_length = VA::size();
301 const unsigned int n_owned = offline_data_->n_locally_owned();
302 const unsigned int n_regular = n_owned / simd_length * simd_length;
304 const auto &sparsity_simd = offline_data_->sparsity_pattern_simd();
309 std::cout <<
" perform time-step with tau = " << tau << std::endl;
310 if (crank_nicolson_extrapolation)
311 std::cout <<
" and extrapolate to t + 2 * tau" << std::endl;
327 std::atomic<bool> restart_needed =
false;
338 std::atomic<bool> correction_needed =
false;
348 Scope scope(computing_timer_,
"time step [P] 1 - update velocities");
352 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
353 using T =
decltype(sentinel);
354 unsigned int stride_size = get_stride_size<T>;
356 const auto view = hyperbolic_system_->template view<dim, T>();
359 for (
unsigned int i = left; i < right; i += stride_size) {
360 const auto U_i = old_U.template get_tensor<T>(i);
361 const auto rho_i = view.density(U_i);
362 const auto M_i = view.momentum(U_i);
363 const auto rho_e_i = view.internal_energy(U_i);
364 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
366 write_entry<T>(density_, rho_i, i);
368 for (
unsigned int d = 0; d < dim; ++d) {
369 write_entry<T>(velocity_.block(d), M_i[d] / rho_i, i);
370 write_entry<T>(velocity_rhs_.block(d), m_i * (M_i[d]), i);
372 write_entry<T>(internal_energy_, rho_e_i / rho_i, i);
377 loop(Number(), n_regular, n_owned);
379 loop(VA(), 0, n_regular);
390 const auto &boundary_map = offline_data_->boundary_map();
392 for (
auto entry : boundary_map) {
394 const auto i = std::get<0>(entry);
398 const auto normal = std::get<1>(entry);
399 const auto id = std::get<4>(entry);
400 const auto position = std::get<5>(entry);
404 Tensor<1, dim, Number> V_i;
405 Tensor<1, dim, Number> RHS_i;
406 for (
unsigned int d = 0; d < dim; ++d) {
407 V_i[d] = velocity_.block(d).local_element(i);
408 RHS_i[d] = velocity_rhs_.block(d).local_element(i);
410 V_i -= 1. * (V_i * normal) * normal;
411 RHS_i -= 1. * (RHS_i * normal) * normal;
412 for (
unsigned int d = 0; d < dim; ++d) {
413 velocity_.block(d).local_element(i) = V_i[d];
414 velocity_rhs_.block(d).local_element(i) = RHS_i[d];
420 for (
unsigned int d = 0; d < dim; ++d) {
421 velocity_.block(d).local_element(i) = Number(0.);
422 velocity_rhs_.block(d).local_element(i) = Number(0.);
428 const auto U_i = initial_values_->initial_state(position, t + tau);
429 const auto view = hyperbolic_system_->template view<dim, Number>();
430 const auto rho_i = view.density(U_i);
431 const auto V_i = view.momentum(U_i) / rho_i;
432 const auto e_i = view.internal_energy(U_i) / rho_i;
434 for (
unsigned int d = 0; d < dim; ++d) {
435 velocity_.block(d).local_element(i) = V_i[d];
436 velocity_rhs_.block(d).local_element(i) = V_i[d];
439 internal_energy_.local_element(i) = e_i;
450 affine_constraints.set_zero(density_);
451 affine_constraints.set_zero(internal_energy_);
452 for (
unsigned int d = 0; d < dim; ++d) {
453 affine_constraints.set_zero(velocity_.block(d));
454 affine_constraints.set_zero(velocity_rhs_.block(d));
460 lumped_mass_matrix, density_, affine_constraints);
467 if (use_gmg_velocity_ && reinitialize_gmg) {
468 MGLevelObject<
typename PreconditionChebyshev<
469 VelocityMatrix<dim, float, Number>,
470 LinearAlgebra::distributed::BlockVector<float>,
471 DiagonalMatrix<dim, float>>::AdditionalData>
472 smoother_data(level_matrix_free_.min_level(),
473 level_matrix_free_.max_level());
475 level_velocity_matrices_.resize(level_matrix_free_.min_level(),
476 level_matrix_free_.max_level());
477 mg_transfer_velocity_.interpolate_to_mg(
478 offline_data_->dof_handler(), level_density_, density_);
480 for (
unsigned int level = level_matrix_free_.min_level();
481 level <= level_matrix_free_.max_level();
483 level_velocity_matrices_[level].initialize(
486 level_matrix_free_[level],
487 level_density_[level],
490 level_velocity_matrices_[level].compute_diagonal(
491 smoother_data[level].preconditioner);
492 if (level == level_matrix_free_.min_level()) {
493 smoother_data[level].degree = numbers::invalid_unsigned_int;
494 smoother_data[level].eig_cg_n_iterations = 500;
495 smoother_data[level].smoothing_range = 1e-3;
497 smoother_data[level].degree = gmg_smoother_degree_;
498 smoother_data[level].eig_cg_n_iterations =
499 gmg_smoother_n_cg_iter_;
500 smoother_data[level].smoothing_range = gmg_smoother_range_vel_;
501 if (gmg_smoother_n_cg_iter_ == 0)
502 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_vel_;
505 mg_smoother_velocity_.initialize(level_velocity_matrices_,
515 Scope scope(computing_timer_,
516 "time step [P] _ - synchronization barriers");
522 *std::min_element(internal_energy_.begin(), internal_energy_.end());
524 e_min_old = Utilities::MPI::min(e_min_old,
525 mpi_ensemble_.ensemble_communicator());
528 constexpr Number eps = std::numeric_limits<Number>::epsilon();
529 e_min_old *= (1. - 1000. * eps);
536 Scope scope(computing_timer_,
"time step [P] 1 - update velocities");
540 VelocityMatrix<dim, Number, Number> velocity_operator;
541 velocity_operator.initialize(
542 *parabolic_system_, *offline_data_, matrix_free_, density_, tau);
544 const auto tolerance_velocity =
545 (tolerance_linfty_norm_ ? velocity_rhs_.linfty_norm()
546 : velocity_rhs_.l2_norm()) *
555 if (!use_gmg_velocity_)
556 throw SolverControl::NoConvergence(0, 0.);
558 using bvt_float = LinearAlgebra::distributed::BlockVector<float>;
560 MGCoarseGridApplySmoother<bvt_float> mg_coarse;
561 mg_coarse.initialize(mg_smoother_velocity_);
563 mg::Matrix<bvt_float> mg_matrix(level_velocity_matrices_);
565 Multigrid<bvt_float> mg(mg_matrix,
567 mg_transfer_velocity_,
568 mg_smoother_velocity_,
569 mg_smoother_velocity_,
570 level_velocity_matrices_.min_level(),
571 level_velocity_matrices_.max_level());
573 const auto &dof_handler = offline_data_->dof_handler();
574 PreconditionMG<dim, bvt_float, MGTransferVelocity<dim, float>>
575 preconditioner(dof_handler, mg, mg_transfer_velocity_);
577 SolverControl solver_control(gmg_max_iter_vel_, tolerance_velocity);
578 SolverCG<BlockVector> solver(solver_control);
580 velocity_operator, velocity_, velocity_rhs_, preconditioner);
583 n_iterations_velocity_ =
584 0.9 * n_iterations_velocity_ + 0.1 * solver_control.last_step();
586 }
catch (SolverControl::NoConvergence &) {
588 SolverControl solver_control(1000, tolerance_velocity);
589 SolverCG<BlockVector> solver(solver_control);
591 velocity_operator, velocity_, velocity_rhs_, diagonal_matrix);
594 n_iterations_velocity_ *= 0.9;
595 n_iterations_velocity_ +=
596 0.1 * (use_gmg_velocity_ ? gmg_max_iter_vel_ : 0) +
597 0.1 * solver_control.last_step();
607 Scope scope(computing_timer_,
608 "time step [P] 2 - update internal energy");
613 matrix_free_.template cell_loop<ScalarVector, BlockVector>(
614 [
this](
const auto &data,
617 const auto cell_range) {
618 FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
620 FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(data);
622 const auto mu = parabolic_system_->mu();
623 const auto lambda = parabolic_system_->lambda();
625 for (
unsigned int cell = cell_range.first;
626 cell < cell_range.second;
628 velocity.reinit(cell);
630 velocity.gather_evaluate(src, EvaluationFlags::gradients);
632 for (
unsigned int q = 0; q < velocity.n_q_points; ++q) {
633 if constexpr (dim == 1) {
635 const auto gradient = velocity.get_gradient(q);
636 auto S = (4. / 3. * mu + lambda) * gradient;
637 energy.submit_value(gradient * S, q);
641 const auto symmetric_gradient =
642 velocity.get_symmetric_gradient(q);
643 const auto divergence = trace(symmetric_gradient);
644 auto S = 2. * mu * symmetric_gradient;
645 for (
unsigned int d = 0; d < dim; ++d)
646 S[d][d] += (lambda - 2. / 3. * mu) * divergence;
647 energy.submit_value(symmetric_gradient * S, q);
650 energy.integrate_scatter(EvaluationFlags::values, dst);
653 internal_energy_rhs_,
657 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
661 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
662 using T =
decltype(sentinel);
663 unsigned int stride_size = get_stride_size<T>;
665 const auto view = hyperbolic_system_->template view<dim, T>();
668 for (
unsigned int i = left; i < right; i += stride_size) {
669 const auto rhs_i = get_entry<T>(internal_energy_rhs_, i);
670 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
671 const auto rho_i = get_entry<T>(density_, i);
672 const auto e_i = get_entry<T>(internal_energy_, i);
674 const auto U_i = old_U.template get_tensor<T>(i);
675 const auto V_i = view.momentum(U_i) / rho_i;
677 dealii::Tensor<1, dim, T> V_i_new;
678 for (
unsigned int d = 0; d < dim; ++d) {
679 V_i_new[d] = get_entry<T>(velocity_.block(d), i);
686 const auto correction =
687 crank_nicolson_extrapolation
689 : Number(0.5) * (V_i - V_i_new).norm_square();
692 const auto result = m_i * rho_i * (e_i + correction) + tau * rhs_i;
693 write_entry<T>(internal_energy_rhs_, result, i);
698 loop(Number(), n_regular, n_owned);
700 loop(VA(), 0, n_regular);
711 const auto &boundary_map = offline_data_->boundary_map();
713 for (
auto entry : boundary_map) {
715 const auto i = std::get<0>(entry);
719 const auto id = std::get<4>(entry);
720 const auto position = std::get<5>(entry);
724 const auto U_i = initial_values_->initial_state(position, t + tau);
725 const auto view = hyperbolic_system_->template view<dim, Number>();
726 const auto rho_i = view.density(U_i);
727 const auto e_i = view.internal_energy(U_i) / rho_i;
728 internal_energy_rhs_.local_element(i) = e_i;
738 affine_constraints.set_zero(internal_energy_);
739 affine_constraints.set_zero(internal_energy_rhs_);
746 if (use_gmg_internal_energy_ && reinitialize_gmg) {
747 MGLevelObject<
typename PreconditionChebyshev<
748 EnergyMatrix<dim, float, Number>,
749 LinearAlgebra::distributed::Vector<float>>::AdditionalData>
750 smoother_data(level_matrix_free_.min_level(),
751 level_matrix_free_.max_level());
753 level_energy_matrices_.resize(level_matrix_free_.min_level(),
754 level_matrix_free_.max_level());
756 for (
unsigned int level = level_matrix_free_.min_level();
757 level <= level_matrix_free_.max_level();
759 level_energy_matrices_[level].initialize(
761 level_matrix_free_[level],
762 level_density_[level],
763 tau * parabolic_system_->cv_inverse_kappa(),
765 level_energy_matrices_[level].compute_diagonal(
766 smoother_data[level].preconditioner);
767 if (level == level_matrix_free_.min_level()) {
768 smoother_data[level].degree = numbers::invalid_unsigned_int;
769 smoother_data[level].eig_cg_n_iterations = 500;
770 smoother_data[level].smoothing_range = 1e-3;
772 smoother_data[level].degree = gmg_smoother_degree_;
773 smoother_data[level].eig_cg_n_iterations =
774 gmg_smoother_n_cg_iter_;
775 smoother_data[level].smoothing_range = gmg_smoother_range_en_;
776 if (gmg_smoother_n_cg_iter_ == 0)
777 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_en_;
780 mg_smoother_energy_.initialize(level_energy_matrices_, smoother_data);
790 Scope scope(computing_timer_,
791 "time step [P] 2 - update internal energy");
795 EnergyMatrix<dim, Number, Number> energy_operator;
796 const auto &kappa = parabolic_system_->cv_inverse_kappa();
797 energy_operator.initialize(
798 *offline_data_, matrix_free_, density_, tau * kappa);
800 const auto tolerance_internal_energy =
801 (tolerance_linfty_norm_ ? internal_energy_rhs_.linfty_norm()
802 : internal_energy_rhs_.l2_norm()) *
806 if (!use_gmg_internal_energy_)
807 throw SolverControl::NoConvergence(0, 0.);
809 using vt_float = LinearAlgebra::distributed::Vector<float>;
810 MGCoarseGridApplySmoother<vt_float> mg_coarse;
811 mg_coarse.initialize(mg_smoother_energy_);
812 mg::Matrix<vt_float> mg_matrix(level_energy_matrices_);
814 Multigrid<vt_float> mg(mg_matrix,
819 level_energy_matrices_.min_level(),
820 level_energy_matrices_.max_level());
822 const auto &dof_handler = offline_data_->dof_handler();
823 PreconditionMG<dim, vt_float, MGTransferEnergy<dim, float>>
824 preconditioner(dof_handler, mg, mg_transfer_energy_);
826 SolverControl solver_control(gmg_max_iter_en_,
827 tolerance_internal_energy);
828 SolverCG<ScalarVector> solver(solver_control);
829 solver.solve(energy_operator,
831 internal_energy_rhs_,
835 n_iterations_internal_energy_ = 0.9 * n_iterations_internal_energy_ +
836 0.1 * solver_control.last_step();
838 }
catch (SolverControl::NoConvergence &) {
840 SolverControl solver_control(1000, tolerance_internal_energy);
841 SolverCG<ScalarVector> solver(solver_control);
842 solver.solve(energy_operator,
844 internal_energy_rhs_,
848 n_iterations_internal_energy_ *= 0.9;
849 n_iterations_internal_energy_ +=
850 0.1 * (use_gmg_internal_energy_ ? gmg_max_iter_en_ : 0) +
851 0.1 * solver_control.last_step();
864 Scope scope(computing_timer_,
"time step [P] 3 - write back vectors");
869 auto loop = [&](
auto sentinel,
unsigned int left,
unsigned int right) {
870 using T =
decltype(sentinel);
871 unsigned int stride_size = get_stride_size<T>;
873 const auto view = hyperbolic_system_->template view<dim, T>();
876 for (
unsigned int i = left; i < right; i += stride_size) {
879 const unsigned int row_length = sparsity_simd.row_length(i);
883 auto U_i = old_U.template get_tensor<T>(i);
884 const auto rho_i = view.density(U_i);
886 Tensor<1, dim, T> m_i_new;
887 for (
unsigned int d = 0; d < dim; ++d) {
888 m_i_new[d] = rho_i * get_entry<T>(velocity_.block(d), i);
891 auto rho_e_i_new = rho_i * get_entry<T>(internal_energy_, i);
898 if (!(T(0.) == std::max(T(0.), rho_i * e_min_old - rho_e_i_new))) {
900 std::cout << std::fixed << std::setprecision(16);
901 const auto e_i_new = rho_e_i_new / rho_i;
902 std::cout <<
"Bounds violation: internal energy (critical)!\n"
903 <<
"\t\te_min_old: " << e_min_old <<
"\n"
904 <<
"\t\te_min_old (delta): "
906 <<
"\t\te_min_new: " << e_i_new <<
"\n"
909 restart_needed =
true;
912 if (crank_nicolson_extrapolation) {
913 m_i_new = Number(2.0) * m_i_new - view.momentum(U_i);
915 Number(2.0) * rho_e_i_new - view.internal_energy(U_i);
923 std::max(T(0.), eps * rho_i * e_min_old - rho_e_i_new))) {
925 std::cout << std::fixed << std::setprecision(16);
926 const auto e_i_new = rho_e_i_new / rho_i;
928 std::cout <<
"Bounds violation: high-order internal energy!"
929 <<
"\t\te_min_new: " << e_i_new <<
"\n"
930 <<
"\t\t-- correction required --" << std::endl;
932 correction_needed =
true;
936 const auto E_i_new = rho_e_i_new + 0.5 * m_i_new * m_i_new / rho_i;
938 for (
unsigned int d = 0; d < dim; ++d)
939 U_i[1 + d] = m_i_new[d];
940 U_i[1 + dim] = E_i_new;
942 new_U.template write_tensor<T>(U_i, i);
947 loop(Number(), n_regular, n_owned);
949 loop(VA(), 0, n_regular);
954 new_U.update_ghost_values();
960 Scope scope(computing_timer_,
961 "time step [H] _ - synchronization barriers");
972 restart_needed.store(Utilities::MPI::logical_or(
973 restart_needed.load(),
974 mpi_ensemble_.synchronization_communicator()));
976 correction_needed.store(Utilities::MPI::logical_or(
977 correction_needed.load(),
978 mpi_ensemble_.synchronization_communicator()));
981 if (correction_needed) {
992 if (restart_needed) {
993 switch (id_violation_strategy) {
1005 template <
typename Description,
int dim,
typename Number>
1007 std::ostream &output)
const
1009 output <<
" [ " << std::setprecision(2) << std::fixed
1010 << n_iterations_velocity_
1011 << (use_gmg_velocity_ ?
" GMG vel -- " :
" CG vel -- ")
1012 << n_iterations_internal_energy_
1013 << (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
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