ryujin 2.1.1 revision d1a5601757449924e68a428cfd892dfe8915810d
parabolic_solver.template.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2023 by the ryujin authors
4//
5
6#pragma once
7
8#include "parabolic_solver.h"
9
10#include <introspection.h>
11#include <openmp.h>
12#include <scope.h>
13#include <simd.h>
14
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>
24
25#include <atomic>
26
27namespace ryujin
28{
29 namespace NavierStokes
30 {
31 using namespace dealii;
32
33 template <typename Description, int dim, typename Number>
35 const MPI_Comm &mpi_communicator,
36 std::map<std::string, dealii::Timer> &computing_timer,
37 const HyperbolicSystem &hyperbolic_system,
38 const ParabolicSystem &parabolic_system,
39 const OfflineData<dim, Number> &offline_data,
40 const InitialValues<Description, dim, Number> &initial_values,
41 const std::string &subsection /*= "ParabolicSolver"*/)
42 : ParameterAcceptor(subsection)
43 , mpi_communicator_(mpi_communicator)
44 , computing_timer_(computing_timer)
45 , hyperbolic_system_(&hyperbolic_system)
46 , parabolic_system_(&parabolic_system)
47 , offline_data_(&offline_data)
48 , initial_values_(&initial_values)
49 , n_restarts_(0)
50 , n_warnings_(0)
51 , n_iterations_velocity_(0.)
52 , n_iterations_internal_energy_(0.)
53 {
54 use_gmg_velocity_ = false;
55 add_parameter("multigrid velocity",
56 use_gmg_velocity_,
57 "Use geometric multigrid for velocity component");
58
59 gmg_max_iter_vel_ = 12;
60 add_parameter("multigrid velocity - max iter",
61 gmg_max_iter_vel_,
62 "Maximal number of CG iterations with GMG smoother");
63
64 gmg_smoother_range_vel_ = 8.;
65 add_parameter("multigrid velocity - chebyshev range",
66 gmg_smoother_range_vel_,
67 "Chebyshev smoother: eigenvalue range parameter");
68
69 gmg_smoother_max_eig_vel_ = 2.0;
70 add_parameter("multigrid velocity - chebyshev max eig",
71 gmg_smoother_max_eig_vel_,
72 "Chebyshev smoother: maximal eigenvalue");
73
74 use_gmg_internal_energy_ = false;
75 add_parameter("multigrid energy",
76 use_gmg_internal_energy_,
77 "Use geometric multigrid for internal energy component");
78
79 gmg_max_iter_en_ = 15;
80 add_parameter("multigrid energy - max iter",
81 gmg_max_iter_en_,
82 "Maximal number of CG iterations with GMG smoother");
83
84 gmg_smoother_range_en_ = 15.;
85 add_parameter("multigrid energy - chebyshev range",
86 gmg_smoother_range_en_,
87 "Chebyshev smoother: eigenvalue range parameter");
88
89 gmg_smoother_max_eig_en_ = 2.0;
90 add_parameter("multigrid energy - chebyshev max eig",
91 gmg_smoother_max_eig_en_,
92 "Chebyshev smoother: maximal eigenvalue");
93
94 gmg_smoother_degree_ = 3;
95 add_parameter("multigrid - chebyshev degree",
96 gmg_smoother_degree_,
97 "Chebyshev smoother: degree");
98
99 gmg_smoother_n_cg_iter_ = 10;
100 add_parameter(
101 "multigrid - chebyshev cg iter",
102 gmg_smoother_n_cg_iter_,
103 "Chebyshev smoother: number of CG iterations to approximate "
104 "eigenvalue");
105
106 gmg_min_level_ = 0;
107 add_parameter(
108 "multigrid - min level",
109 gmg_min_level_,
110 "Minimal mesh level to be visited in the geometric multigrid "
111 "cycle where the coarse grid solver (Chebyshev) is called");
112
113 tolerance_ = Number(1.0e-12);
114 add_parameter("tolerance", tolerance_, "Tolerance for linear solvers");
115
116 tolerance_linfty_norm_ = false;
117 add_parameter("tolerance linfty norm",
118 tolerance_linfty_norm_,
119 "Use the l_infty norm instead of the l_2 norm for the "
120 "stopping criterion");
121 }
122
123
124 template <typename Description, int dim, typename Number>
126 {
127#ifdef DEBUG_OUTPUT
128 std::cout << "ParabolicSolver<dim, Number>::prepare()" << std::endl;
129#endif
130
131 const auto &discretization = offline_data_->discretization();
132 AssertThrow(discretization.ansatz() == Ansatz::cg_q1,
133 dealii::ExcMessage("The NavierStokes module currently only "
134 "supports cG Q1 finite elements."));
135
136 /* Initialize vectors: */
137
138 typename MatrixFree<dim, Number>::AdditionalData additional_data;
139 additional_data.tasks_parallel_scheme =
140 MatrixFree<dim, Number>::AdditionalData::none;
141
142 matrix_free_.reinit(discretization.mapping(),
143 offline_data_->dof_handler(),
144 offline_data_->affine_constraints(),
145 discretization.quadrature_1d(),
146 additional_data);
147
148 const auto &scalar_partitioner =
149 matrix_free_.get_dof_info(0).vector_partitioner;
150
151 velocity_.reinit(dim);
152 velocity_rhs_.reinit(dim);
153 for (unsigned int i = 0; i < dim; ++i) {
154 velocity_.block(i).reinit(scalar_partitioner);
155 velocity_rhs_.block(i).reinit(scalar_partitioner);
156 }
157
158 internal_energy_.reinit(scalar_partitioner);
159 internal_energy_rhs_.reinit(scalar_partitioner);
160
161 density_.reinit(scalar_partitioner);
162
163 /* Initialize multigrid: */
164
165 if (!use_gmg_velocity_ && !use_gmg_internal_energy_)
166 return;
167
168 const unsigned int n_levels =
169 offline_data_->dof_handler().get_triangulation().n_global_levels();
170 const unsigned int min_level = std::min(gmg_min_level_, n_levels - 1);
171 MGLevelObject<IndexSet> relevant_sets(0, n_levels - 1);
172 for (unsigned int level = 0; level < n_levels; ++level)
173 dealii::DoFTools::extract_locally_relevant_level_dofs(
174 offline_data_->dof_handler(), level, relevant_sets[level]);
175 mg_constrained_dofs_.initialize(offline_data_->dof_handler(),
176 relevant_sets);
177 std::set<types::boundary_id> boundary_ids;
178 boundary_ids.insert(Boundary::dirichlet);
179 boundary_ids.insert(Boundary::no_slip);
180 mg_constrained_dofs_.make_zero_boundary_constraints(
181 offline_data_->dof_handler(), boundary_ids);
182
183 typename MatrixFree<dim, float>::AdditionalData additional_data_level;
184 additional_data_level.tasks_parallel_scheme =
185 MatrixFree<dim, float>::AdditionalData::none;
186
187 level_matrix_free_.resize(min_level, n_levels - 1);
188 level_density_.resize(min_level, n_levels - 1);
189 for (unsigned int level = min_level; level < n_levels; ++level) {
190 additional_data_level.mg_level = level;
191 AffineConstraints<double> constraints(relevant_sets[level]);
192 // constraints.add_lines(mg_constrained_dofs_.get_boundary_indices(level));
193 // constraints.merge(mg_constrained_dofs_.get_level_constraints(level));
194 constraints.close();
195 level_matrix_free_[level].reinit(discretization.mapping(),
196 offline_data_->dof_handler(),
197 constraints,
198 discretization.quadrature_1d(),
199 additional_data_level);
200 level_matrix_free_[level].initialize_dof_vector(level_density_[level]);
201 }
202
203 mg_transfer_velocity_.build(offline_data_->dof_handler(),
204 mg_constrained_dofs_,
205 level_matrix_free_);
206 mg_transfer_energy_.build(offline_data_->dof_handler(),
207 level_matrix_free_);
208 }
209
210
211 template <typename Description, int dim, typename Number>
213 const StateVector &old_state_vector,
214 const Number t,
215 StateVector &new_state_vector,
216 Number tau,
217 const IDViolationStrategy id_violation_strategy,
218 const bool reinitialize_gmg) const
219 {
220#ifdef DEBUG_OUTPUT
221 std::cout << "ParabolicSolver<dim, Number>::step()" << std::endl;
222#endif
223
224 const auto &old_U = std::get<0>(old_state_vector);
225 auto &new_U = std::get<0>(new_state_vector);
226
228
229 using VA = VectorizedArray<Number>;
230
231 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
232 const auto &affine_constraints = offline_data_->affine_constraints();
233
234 /* Index ranges for the iteration over the sparsity pattern : */
235
236 constexpr auto simd_length = VA::size();
237 const unsigned int n_owned = offline_data_->n_locally_owned();
238 const unsigned int n_regular = n_owned / simd_length * simd_length;
239
240 DiagonalMatrix<dim, Number> diagonal_matrix;
241
242#ifdef DEBUG_OUTPUT
243 std::cout << " perform time-step with tau = " << tau << std::endl;
244#endif
245
246 /* A boolean signalling that a restart is necessary: */
247 std::atomic<bool> restart_needed = false;
248
249 /*
250 * Step 1:
251 *
252 * Build right hand side for the velocity update.
253 * Also initialize solution vectors for internal energy and velocity
254 * update.
255 */
256 {
257 Scope scope(computing_timer_, "time step [P] 1 - update velocities");
259 LIKWID_MARKER_START("time_step_parabolic_1");
260
261 auto loop = [&](auto sentinel, unsigned int left, unsigned int right) {
262 using T = decltype(sentinel);
263 unsigned int stride_size = get_stride_size<T>;
264
265 const auto view = hyperbolic_system_->template view<dim, T>();
266
268 for (unsigned int i = left; i < right; i += stride_size) {
269 const auto U_i = old_U.template get_tensor<T>(i);
270 const auto rho_i = view.density(U_i);
271 const auto M_i = view.momentum(U_i);
272 const auto rho_e_i = view.internal_energy(U_i);
273 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
274
275 write_entry<T>(density_, rho_i, i);
276 /* (5.4a) */
277 for (unsigned int d = 0; d < dim; ++d) {
278 write_entry<T>(velocity_.block(d), M_i[d] / rho_i, i);
279 write_entry<T>(velocity_rhs_.block(d), m_i * (M_i[d]), i);
280 }
281 write_entry<T>(internal_energy_, rho_e_i / rho_i, i);
282 }
283 };
284
285 /* Parallel non-vectorized loop: */
286 loop(Number(), n_regular, n_owned);
287 /* Parallel vectorized SIMD loop: */
288 loop(VA(), 0, n_regular);
289
291
292 /*
293 * Set up "strongly enforced" boundary conditions that are not stored
294 * in the AffineConstraints map. In this case we enforce boundary
295 * values by imposing them strongly in the iteration by setting the
296 * initial vector and the right hand side to the right value:
297 */
298
299 const auto &boundary_map = offline_data_->boundary_map();
300
301 for (auto entry : boundary_map) {
302 const auto i = entry.first;
303 if (i >= n_owned)
304 continue;
305
306 const auto normal = std::get<0>(entry.second);
307 const auto id = std::get<3>(entry.second);
308 const auto position = std::get<4>(entry.second);
309
310 if (id == Boundary::slip) {
311 /* Remove normal component of velocity: */
312 Tensor<1, dim, Number> V_i;
313 Tensor<1, dim, Number> RHS_i;
314 for (unsigned int d = 0; d < dim; ++d) {
315 V_i[d] = velocity_.block(d).local_element(i);
316 RHS_i[d] = velocity_rhs_.block(d).local_element(i);
317 }
318 V_i -= 1. * (V_i * normal) * normal;
319 RHS_i -= 1. * (RHS_i * normal) * normal;
320 for (unsigned int d = 0; d < dim; ++d) {
321 velocity_.block(d).local_element(i) = V_i[d];
322 velocity_rhs_.block(d).local_element(i) = RHS_i[d];
323 }
324
325 } else if (id == Boundary::no_slip) {
326
327 /* Set velocity to zero: */
328 for (unsigned int d = 0; d < dim; ++d) {
329 velocity_.block(d).local_element(i) = Number(0.);
330 velocity_rhs_.block(d).local_element(i) = Number(0.);
331 }
332
333 } else if (id == Boundary::dirichlet) {
334
335 /* Prescribe velocity: */
336 const auto U_i = initial_values_->initial_state(position, t + tau);
337 const auto view = hyperbolic_system_->template view<dim, Number>();
338 const auto rho_i = view.density(U_i);
339 const auto V_i = view.momentum(U_i) / rho_i;
340 const auto e_i = view.internal_energy(U_i) / rho_i;
341
342 for (unsigned int d = 0; d < dim; ++d) {
343 velocity_.block(d).local_element(i) = V_i[d];
344 velocity_rhs_.block(d).local_element(i) = V_i[d];
345 }
346
347 internal_energy_.local_element(i) = e_i;
348 }
349 }
350
351 /*
352 * Zero out constrained degrees of freedom due to periodic boundary
353 * conditions. These boundary conditions are enforced by modifying
354 * the stencil - consequently we have to remove constrained dofs from
355 * the linear system.
356 */
357
358 affine_constraints.set_zero(density_);
359 affine_constraints.set_zero(internal_energy_);
360 for (unsigned int d = 0; d < dim; ++d) {
361 affine_constraints.set_zero(velocity_.block(d));
362 affine_constraints.set_zero(velocity_rhs_.block(d));
363 }
364
365 /* Prepare preconditioner: */
366
367 diagonal_matrix.reinit(
368 lumped_mass_matrix, density_, affine_constraints);
369
370 /*
371 * Update MG matrices all 4 time steps; this is a balance because more
372 * refreshes will render the approximation better, at some additional
373 * cost.
374 */
375 if (use_gmg_velocity_ && reinitialize_gmg) {
376 MGLevelObject<typename PreconditionChebyshev<
378 LinearAlgebra::distributed::BlockVector<float>,
379 DiagonalMatrix<dim, float>>::AdditionalData>
380 smoother_data(level_matrix_free_.min_level(),
381 level_matrix_free_.max_level());
382
383 level_velocity_matrices_.resize(level_matrix_free_.min_level(),
384 level_matrix_free_.max_level());
385 mg_transfer_velocity_.interpolate_to_mg(
386 offline_data_->dof_handler(), level_density_, density_);
387
388 for (unsigned int level = level_matrix_free_.min_level();
389 level <= level_matrix_free_.max_level();
390 ++level) {
391 level_velocity_matrices_[level].initialize(
392 *parabolic_system_,
393 *offline_data_,
394 level_matrix_free_[level],
395 level_density_[level],
396 tau,
397 level);
398 level_velocity_matrices_[level].compute_diagonal(
399 smoother_data[level].preconditioner);
400 if (level == level_matrix_free_.min_level()) {
401 smoother_data[level].degree = numbers::invalid_unsigned_int;
402 smoother_data[level].eig_cg_n_iterations = 500;
403 smoother_data[level].smoothing_range = 1e-3;
404 } else {
405 smoother_data[level].degree = gmg_smoother_degree_;
406 smoother_data[level].eig_cg_n_iterations =
407 gmg_smoother_n_cg_iter_;
408 smoother_data[level].smoothing_range = gmg_smoother_range_vel_;
409 if (gmg_smoother_n_cg_iter_ == 0)
410 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_vel_;
411 }
412 }
413 mg_smoother_velocity_.initialize(level_velocity_matrices_,
414 smoother_data);
415 }
416
417 LIKWID_MARKER_STOP("time_step_parabolic_1");
418 }
419
420 /* Compute the global minimum of the internal energy: */
421
422 // .begin() and .end() denote the locally owned index range:
423 auto e_min_old =
424 *std::min_element(internal_energy_.begin(), internal_energy_.end());
425 e_min_old = Utilities::MPI::min(e_min_old, mpi_communicator_);
426
427 // FIXME: create a meaningful relaxation based on global mesh size min.
428 constexpr Number eps = std::numeric_limits<Number>::epsilon();
429 e_min_old *= (1. - 1000. * eps);
430
431 /*
432 * Step 1: Solve velocity update:
433 */
434 {
435 Scope scope(computing_timer_, "time step [P] 1 - update velocities");
436
437 LIKWID_MARKER_START("time_step_parabolic_1");
438
439 VelocityMatrix<dim, Number, Number> velocity_operator;
440 velocity_operator.initialize(
441 *parabolic_system_, *offline_data_, matrix_free_, density_, tau);
442
443 const auto tolerance_velocity =
444 (tolerance_linfty_norm_ ? velocity_rhs_.linfty_norm()
445 : velocity_rhs_.l2_norm()) *
446 tolerance_;
447
448 /*
449 * Multigrid might lack robustness for some cases, so in case it takes
450 * too many iterations we better switch to the more robust plain
451 * conjugate gradient method.
452 */
453 try {
454 if (!use_gmg_velocity_)
455 throw SolverControl::NoConvergence(0, 0.);
456
457 using bvt_float = LinearAlgebra::distributed::BlockVector<float>;
458
459 MGCoarseGridApplySmoother<bvt_float> mg_coarse;
460 mg_coarse.initialize(mg_smoother_velocity_);
461
462 mg::Matrix<bvt_float> mg_matrix(level_velocity_matrices_);
463
464 Multigrid<bvt_float> mg(mg_matrix,
465 mg_coarse,
466 mg_transfer_velocity_,
467 mg_smoother_velocity_,
468 mg_smoother_velocity_,
469 level_velocity_matrices_.min_level(),
470 level_velocity_matrices_.max_level());
471
472 const auto &dof_handler = offline_data_->dof_handler();
473 PreconditionMG<dim, bvt_float, MGTransferVelocity<dim, float>>
474 preconditioner(dof_handler, mg, mg_transfer_velocity_);
475
476 SolverControl solver_control(gmg_max_iter_vel_, tolerance_velocity);
477 SolverCG<BlockVector> solver(solver_control);
478 solver.solve(
479 velocity_operator, velocity_, velocity_rhs_, preconditioner);
480
481 /* update exponential moving average */
482 n_iterations_velocity_ =
483 0.9 * n_iterations_velocity_ + 0.1 * solver_control.last_step();
484
485 } catch (SolverControl::NoConvergence &) {
486
487 SolverControl solver_control(1000, tolerance_velocity);
488 SolverCG<BlockVector> solver(solver_control);
489 solver.solve(
490 velocity_operator, velocity_, velocity_rhs_, diagonal_matrix);
491
492 /* update exponential moving average, counting also GMG iterations */
493 n_iterations_velocity_ *= 0.9;
494 n_iterations_velocity_ +=
495 0.1 * (use_gmg_velocity_ ? gmg_max_iter_vel_ : 0) +
496 0.1 * solver_control.last_step();
497 }
498
499 LIKWID_MARKER_STOP("time_step_parabolic_1");
500 }
501
502 /*
503 * Step 2: Build internal energy right hand side:
504 */
505 {
506 Scope scope(computing_timer_,
507 "time step [P] 2 - update internal energy");
508
509 LIKWID_MARKER_START("time_step_parabolic_2");
510
511 /* Compute m_i K_i^{n+1/2}: (5.5) */
512 matrix_free_.template cell_loop<ScalarVector, BlockVector>(
513 [this](const auto &data,
514 auto &dst,
515 const auto &src,
516 const auto cell_range) {
517 FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
518 data);
519 FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(data);
520
521 const auto mu = parabolic_system_->mu();
522 const auto lambda = parabolic_system_->lambda();
523
524 for (unsigned int cell = cell_range.first;
525 cell < cell_range.second;
526 ++cell) {
527 velocity.reinit(cell);
528 energy.reinit(cell);
529 velocity.gather_evaluate(src, EvaluationFlags::gradients);
530
531 for (unsigned int q = 0; q < velocity.n_q_points; ++q) {
532 if constexpr (dim == 1) {
533 /* Workaround: no symmetric gradient for dim == 1: */
534 const auto gradient = velocity.get_gradient(q);
535 auto S = (4. / 3. * mu + lambda) * gradient;
536 energy.submit_value(gradient * S, q);
537
538 } else {
539
540 const auto symmetric_gradient =
541 velocity.get_symmetric_gradient(q);
542 const auto divergence = trace(symmetric_gradient);
543 auto S = 2. * mu * symmetric_gradient;
544 for (unsigned int d = 0; d < dim; ++d)
545 S[d][d] += (lambda - 2. / 3. * mu) * divergence;
546 energy.submit_value(symmetric_gradient * S, q);
547 }
548 }
549 energy.integrate_scatter(EvaluationFlags::values, dst);
550 }
551 },
552 internal_energy_rhs_,
553 velocity_,
554 /* zero destination */ true);
555
556 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
557
559
560 auto loop = [&](auto sentinel, unsigned int left, unsigned int right) {
561 using T = decltype(sentinel);
562 unsigned int stride_size = get_stride_size<T>;
563
564 const auto view = hyperbolic_system_->template view<dim, T>();
565
567 for (unsigned int i = left; i < right; i += stride_size) {
568 const auto rhs_i = get_entry<T>(internal_energy_rhs_, i);
569 const auto m_i = get_entry<T>(lumped_mass_matrix, i);
570 const auto rho_i = get_entry<T>(density_, i);
571 const auto e_i = get_entry<T>(internal_energy_, i);
572
573 const auto U_i = old_U.template get_tensor<T>(i);
574 const auto V_i = view.momentum(U_i) / rho_i;
575
576 dealii::Tensor<1, dim, T> V_i_new;
577 for (unsigned int d = 0; d < dim; ++d) {
578 V_i_new[d] = get_entry<T>(velocity_.block(d), i);
579 }
580
581 /*
582 * For backward Euler we have to add this algebraic correction
583 * to ensure conservation of total energy.
584 */
585 const auto correction = Number(0.5) * (V_i - V_i_new).norm_square();
586
587 /* rhs_i contains already m_i K_i^{n+1/2} */
588 const auto result = m_i * rho_i * (e_i + correction) + tau * rhs_i;
589 write_entry<T>(internal_energy_rhs_, result, i);
590 }
591 };
592
593 /* Parallel non-vectorized loop: */
594 loop(Number(), n_regular, n_owned);
595 /* Parallel vectorized SIMD loop: */
596 loop(VA(), 0, n_regular);
597
599
600 /*
601 * Set up "strongly enforced" boundary conditions that are not stored
602 * in the AffineConstraints map: We enforce Neumann conditions (i.e.,
603 * insulating boundary conditions) everywhere except for Dirichlet
604 * boundaries where we have to enforce prescribed conditions:
605 */
606
607 const auto &boundary_map = offline_data_->boundary_map();
608
609 for (auto entry : boundary_map) {
610 const auto i = entry.first;
611 if (i >= n_owned)
612 continue;
613
614 const auto id = std::get<3>(entry.second);
615 const auto position = std::get<4>(entry.second);
616
617 if (id == Boundary::dirichlet) {
618 /* Prescribe internal energy: */
619 const auto U_i = initial_values_->initial_state(position, t + tau);
620 const auto view = hyperbolic_system_->template view<dim, Number>();
621 const auto rho_i = view.density(U_i);
622 const auto e_i = view.internal_energy(U_i) / rho_i;
623 internal_energy_rhs_.local_element(i) = e_i;
624 }
625 }
626
627 /*
628 * Zero out constrained degrees of freedom due to periodic boundary
629 * conditions. These boundary conditions are enforced by modifying
630 * the stencil - consequently we have to remove constrained dofs from
631 * the linear system.
632 */
633 affine_constraints.set_zero(internal_energy_rhs_);
634
635 /*
636 * Update MG matrices all 4 time steps; this is a balance because more
637 * refreshes will render the approximation better, at some additional
638 * cost.
639 */
640 if (use_gmg_internal_energy_ && reinitialize_gmg) {
641 MGLevelObject<typename PreconditionChebyshev<
643 LinearAlgebra::distributed::Vector<float>>::AdditionalData>
644 smoother_data(level_matrix_free_.min_level(),
645 level_matrix_free_.max_level());
646
647 level_energy_matrices_.resize(level_matrix_free_.min_level(),
648 level_matrix_free_.max_level());
649
650 for (unsigned int level = level_matrix_free_.min_level();
651 level <= level_matrix_free_.max_level();
652 ++level) {
653 level_energy_matrices_[level].initialize(
654 *offline_data_,
655 level_matrix_free_[level],
656 level_density_[level],
657 tau * parabolic_system_->cv_inverse_kappa(),
658 level);
659 level_energy_matrices_[level].compute_diagonal(
660 smoother_data[level].preconditioner);
661 if (level == level_matrix_free_.min_level()) {
662 smoother_data[level].degree = numbers::invalid_unsigned_int;
663 smoother_data[level].eig_cg_n_iterations = 500;
664 smoother_data[level].smoothing_range = 1e-3;
665 } else {
666 smoother_data[level].degree = gmg_smoother_degree_;
667 smoother_data[level].eig_cg_n_iterations =
668 gmg_smoother_n_cg_iter_;
669 smoother_data[level].smoothing_range = gmg_smoother_range_en_;
670 if (gmg_smoother_n_cg_iter_ == 0)
671 smoother_data[level].max_eigenvalue = gmg_smoother_max_eig_en_;
672 }
673 }
674 mg_smoother_energy_.initialize(level_energy_matrices_, smoother_data);
675 }
676
677 LIKWID_MARKER_STOP("time_step_parabolic_2");
678 }
679
680 /*
681 * Step 2: Solve internal energy update:
682 */
683 {
684 Scope scope(computing_timer_,
685 "time step [P] 2 - update internal energy");
686
687 LIKWID_MARKER_START("time_step_parabolic_2");
688
689 EnergyMatrix<dim, Number, Number> energy_operator;
690 const auto &kappa = parabolic_system_->cv_inverse_kappa();
691 energy_operator.initialize(
692 *offline_data_, matrix_free_, density_, tau * kappa);
693
694 const auto tolerance_internal_energy =
695 (tolerance_linfty_norm_ ? internal_energy_rhs_.linfty_norm()
696 : internal_energy_rhs_.l2_norm()) *
697 tolerance_;
698
699 try {
700 if (!use_gmg_internal_energy_)
701 throw SolverControl::NoConvergence(0, 0.);
702
703 using vt_float = LinearAlgebra::distributed::Vector<float>;
704 MGCoarseGridApplySmoother<vt_float> mg_coarse;
705 mg_coarse.initialize(mg_smoother_energy_);
706 mg::Matrix<vt_float> mg_matrix(level_energy_matrices_);
707
708 Multigrid<vt_float> mg(mg_matrix,
709 mg_coarse,
710 mg_transfer_energy_,
711 mg_smoother_energy_,
712 mg_smoother_energy_,
713 level_energy_matrices_.min_level(),
714 level_energy_matrices_.max_level());
715
716 const auto &dof_handler = offline_data_->dof_handler();
717 PreconditionMG<dim, vt_float, MGTransferEnergy<dim, float>>
718 preconditioner(dof_handler, mg, mg_transfer_energy_);
719
720 SolverControl solver_control(gmg_max_iter_en_,
721 tolerance_internal_energy);
722 SolverCG<ScalarVector> solver(solver_control);
723 solver.solve(energy_operator,
724 internal_energy_,
725 internal_energy_rhs_,
726 preconditioner);
727
728 /* update exponential moving average */
729 n_iterations_internal_energy_ = 0.9 * n_iterations_internal_energy_ +
730 0.1 * solver_control.last_step();
731
732 } catch (SolverControl::NoConvergence &) {
733
734 SolverControl solver_control(1000, tolerance_internal_energy);
735 SolverCG<ScalarVector> solver(solver_control);
736 solver.solve(energy_operator,
737 internal_energy_,
738 internal_energy_rhs_,
739 diagonal_matrix);
740
741 /* update exponential moving average, counting also GMG iterations */
742 n_iterations_internal_energy_ *= 0.9;
743 n_iterations_internal_energy_ +=
744 0.1 * (use_gmg_internal_energy_ ? gmg_max_iter_en_ : 0) +
745 0.1 * solver_control.last_step();
746 }
747
748 /*
749 * Check for local minimum principle on internal energy:
750 */
751 {
752 // .begin() and .end() denote the locally owned index range:
753 auto e_min_new = *std::min_element(internal_energy_.begin(),
754 internal_energy_.end());
755 e_min_new = Utilities::MPI::min(e_min_new, mpi_communicator_);
756
757 if (e_min_new < e_min_old) {
758#ifdef DEBUG_OUTPUT
759 std::cout << std::fixed << std::setprecision(16);
760 std::cout << "Bounds violation: internal energy (critical)!\n"
761 << "\t\te_min_old: " << e_min_old << "\n"
762 << "\t\te_min_old (delta): "
763 << negative_part(e_min_new - e_min_old) << "\n"
764 << "\t\te_min_new: " << e_min_new << "\n"
765 << std::endl;
766#endif
767 restart_needed = true;
768 }
769 }
770
771 LIKWID_MARKER_STOP("time_step_parabolic_2");
772 }
773
774 /*
775 * Step 3: Copy vectors
776 *
777 * FIXME: Memory access is suboptimal...
778 */
779 {
780 Scope scope(computing_timer_, "time step [P] 3 - write back vectors");
781
783 LIKWID_MARKER_START("time_step_parabolic_3");
784
785 auto loop = [&](auto sentinel, unsigned int left, unsigned int right) {
786 using T = decltype(sentinel);
787 unsigned int stride_size = get_stride_size<T>;
788
789 const auto view = hyperbolic_system_->template view<dim, T>();
790
792 for (unsigned int i = left; i < right; i += stride_size) {
793 auto U_i = old_U.template get_tensor<T>(i);
794 const auto rho_i = view.density(U_i);
795
796 Tensor<1, dim, T> m_i_new;
797 for (unsigned int d = 0; d < dim; ++d) {
798 m_i_new[d] = rho_i * get_entry<T>(velocity_.block(d), i);
799 }
800
801 const auto rho_e_i_new = rho_i * get_entry<T>(internal_energy_, i);
802
803 const auto E_i_new = rho_e_i_new + 0.5 * m_i_new * m_i_new / rho_i;
804
805 for (unsigned int d = 0; d < dim; ++d)
806 U_i[1 + d] = m_i_new[d];
807 U_i[1 + dim] = E_i_new;
808
809 new_U.template write_tensor<T>(U_i, i);
810 }
811 };
812
813 /* Parallel non-vectorized loop: */
814 loop(Number(), n_regular, n_owned);
815 /* Parallel vectorized SIMD loop: */
816 loop(VA(), 0, n_regular);
817
818 LIKWID_MARKER_STOP("time_step_parabolic_3");
820
821 new_U.update_ghost_values();
822 }
823
825
826 if (restart_needed) {
827 switch (id_violation_strategy) {
829 n_warnings_++;
830 break;
832 n_restarts_++;
833 throw Restart();
834 }
835 }
836 }
837
838
839 template <typename Description, int dim, typename Number>
841 std::ostream &output) const
842 {
843 output << " [ " << std::setprecision(2) << std::fixed
844 << n_iterations_velocity_
845 << (use_gmg_velocity_ ? " GMG vel -- " : " CG vel -- ")
846 << n_iterations_internal_energy_
847 << (use_gmg_internal_energy_ ? " GMG int ]" : " CG int ]")
848 << std::endl;
849 }
850
851 } // namespace NavierStokes
852} /* namespace ryujin */
void reinit(const vector_type &lumped_mass_matrix, const vector_type &density, const dealii::AffineConstraints< Number > &affine_constraints)
void initialize(const OfflineData< dim, Number2 > &offline_data, const dealii::MatrixFree< dim, Number > &matrix_free, const dealii::LinearAlgebra::distributed::Vector< Number > &density, const Number time_factor, const unsigned int level=dealii::numbers::invalid_unsigned_int)
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 MPI_Comm &mpi_communicator, std::map< std::string, dealii::Timer > &computing_timer, const HyperbolicSystem &hyperbolic_system, const ParabolicSystem &parabolic_system, const OfflineData< dim, Number > &offline_data, const InitialValues< Description, dim, Number > &initial_values, const std::string &subsection="ParabolicSolver")
void initialize(const ParabolicSystem &parabolic_system, const OfflineData< dim, Number2 > &offline_data, const dealii::MatrixFree< dim, Number > &matrix_free, const dealii::LinearAlgebra::distributed::Vector< Number > &density, const Number theta_x_tau, const unsigned int level=dealii::numbers::invalid_unsigned_int)
#define RYUJIN_PARALLEL_REGION_BEGIN
Definition: openmp.h:54
#define RYUJIN_OMP_FOR
Definition: openmp.h:70
#define RYUJIN_PARALLEL_REGION_END
Definition: openmp.h:63
DEAL_II_ALWAYS_INLINE Number negative_part(const Number number)
Definition: simd.h:124
#define LIKWID_MARKER_START(opt)
Definition: introspection.h:68
#define CALLGRIND_START_INSTRUMENTATION
Definition: introspection.h:28
#define LIKWID_MARKER_STOP(opt)
Definition: introspection.h:73
#define CALLGRIND_STOP_INSTRUMENTATION
Definition: introspection.h:35