ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
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 <instrumentation.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 MPIEnsemble &mpi_ensemble,
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_ensemble_(mpi_ensemble)
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_corrections_(0)
51 , n_warnings_(0)
52 , n_iterations_velocity_(0.)
53 , n_iterations_internal_energy_(0.)
54 {
55 use_gmg_velocity_ = false;
56 add_parameter("multigrid velocity",
57 use_gmg_velocity_,
58 "Use geometric multigrid for velocity component");
59
60 gmg_max_iter_vel_ = 12;
61 add_parameter("multigrid velocity - max iter",
62 gmg_max_iter_vel_,
63 "Maximal number of CG iterations with GMG smoother");
64
65 gmg_smoother_range_vel_ = 8.;
66 add_parameter("multigrid velocity - chebyshev range",
67 gmg_smoother_range_vel_,
68 "Chebyshev smoother: eigenvalue range parameter");
69
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");
74
75 use_gmg_internal_energy_ = false;
76 add_parameter("multigrid energy",
77 use_gmg_internal_energy_,
78 "Use geometric multigrid for internal energy component");
79
80 gmg_max_iter_en_ = 15;
81 add_parameter("multigrid energy - max iter",
82 gmg_max_iter_en_,
83 "Maximal number of CG iterations with GMG smoother");
84
85 gmg_smoother_range_en_ = 15.;
86 add_parameter("multigrid energy - chebyshev range",
87 gmg_smoother_range_en_,
88 "Chebyshev smoother: eigenvalue range parameter");
89
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");
94
95 gmg_smoother_degree_ = 3;
96 add_parameter("multigrid - chebyshev degree",
97 gmg_smoother_degree_,
98 "Chebyshev smoother: degree");
99
100 gmg_smoother_n_cg_iter_ = 10;
101 add_parameter(
102 "multigrid - chebyshev cg iter",
103 gmg_smoother_n_cg_iter_,
104 "Chebyshev smoother: number of CG iterations to approximate "
105 "eigenvalue");
106
107 gmg_min_level_ = 0;
108 add_parameter(
109 "multigrid - min level",
110 gmg_min_level_,
111 "Minimal mesh level to be visited in the geometric multigrid "
112 "cycle where the coarse grid solver (Chebyshev) is called");
113
114 tolerance_ = Number(1.0e-12);
115 add_parameter("tolerance", tolerance_, "Tolerance for linear solvers");
116
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");
122 }
123
124
125 template <typename Description, int dim, typename Number>
127 {
128#ifdef DEBUG_OUTPUT
129 std::cout << "ParabolicSolver<dim, Number>::prepare()" << std::endl;
130#endif
131
132 const auto &discretization = offline_data_->discretization();
133 AssertThrow(discretization.ansatz() == Ansatz::cg_q1,
134 dealii::ExcMessage("The Navier-Stokes module currently only "
135 "supports cG Q1 finite elements."));
136
137 AssertThrow(!offline_data_->dof_handler().has_hp_capabilities(),
138 dealii::ExcMessage(
139 "The Navier-Stokes module currently does not support "
140 "DofHandlers set up with hp capabilities."));
141
142 /* Initialize vectors: */
143
144 typename MatrixFree<dim, Number>::AdditionalData additional_data;
145 additional_data.tasks_parallel_scheme =
146 MatrixFree<dim, Number>::AdditionalData::none;
147
148 matrix_free_.reinit(discretization.mapping(),
149 offline_data_->dof_handler(),
150 offline_data_->affine_constraints(),
151 discretization.quadrature_1d(),
152 additional_data);
153
154 const auto &scalar_partitioner =
155 matrix_free_.get_dof_info(0).vector_partitioner;
156
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);
162 }
163
164 internal_energy_.reinit(scalar_partitioner);
165 internal_energy_rhs_.reinit(scalar_partitioner);
166
167 density_.reinit(scalar_partitioner);
168
169 /* Initialize multigrid: */
170
171 if (!use_gmg_velocity_ && !use_gmg_internal_energy_)
172 return;
173
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);
183#else
184 dealii::DoFTools::extract_locally_relevant_level_dofs(
185 offline_data_->dof_handler(), level, relevant_sets[level]);
186#endif
187 mg_constrained_dofs_.initialize(offline_data_->dof_handler(),
188 relevant_sets);
189 std::set<types::boundary_id> boundary_ids;
190 boundary_ids.insert(Boundary::dirichlet);
191 boundary_ids.insert(Boundary::no_slip);
192 mg_constrained_dofs_.make_zero_boundary_constraints(
193 offline_data_->dof_handler(), boundary_ids);
194
195 typename MatrixFree<dim, float>::AdditionalData additional_data_level;
196 additional_data_level.tasks_parallel_scheme =
197 MatrixFree<dim, float>::AdditionalData::none;
198
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]);
206#else
207 AffineConstraints<double> constraints(relevant_sets[level]);
208#endif
209 // constraints.add_lines(mg_constrained_dofs_.get_boundary_indices(level));
210 // constraints.merge(mg_constrained_dofs_.get_level_constraints(level));
211 constraints.close();
212 level_matrix_free_[level].reinit(discretization.mapping(),
213 offline_data_->dof_handler(),
214 constraints,
215 discretization.quadrature_1d(),
216 additional_data_level);
217 level_matrix_free_[level].initialize_dof_vector(level_density_[level]);
218 }
219
220 mg_transfer_velocity_.build(offline_data_->dof_handler(),
221 mg_constrained_dofs_,
222 level_matrix_free_);
223 mg_transfer_energy_.build(offline_data_->dof_handler(),
224 level_matrix_free_);
225 }
226
227 template <typename Description, int dim, typename Number>
229 StateVector & /*state_vector*/, Number /*t*/) const
230 {
235 }
236
237
238 template <typename Description, int dim, typename Number>
240 const StateVector &old_state_vector,
241 const Number t,
242 StateVector &new_state_vector,
243 Number tau,
244 const IDViolationStrategy id_violation_strategy,
245 const bool reinitialize_gmg) const
246 {
247 /* Backward Euler step to half time step, and extrapolate: */
248
249 step(old_state_vector,
250 t,
251 new_state_vector,
252 tau,
253 id_violation_strategy,
254 reinitialize_gmg,
255 /*crank_nicolson_extrapolation = */ false);
256 }
257
258
259 template <typename Description, int dim, typename Number>
261 const StateVector &old_state_vector,
262 const Number t,
263 StateVector &new_state_vector,
264 Number tau,
265 const IDViolationStrategy id_violation_strategy,
266 const bool reinitialize_gmg) const
267 {
268 try {
269 step(old_state_vector,
270 t,
271 new_state_vector,
272 tau / Number(2.),
273 id_violation_strategy,
274 reinitialize_gmg,
275 /*crank_nicolson_extrapolation = */ true);
276
277 } catch (Correction) {
278
279 /*
280 * Under very rare circumstances we might fail to perform a Crank
281 * Nicolson step because the extrapolation step produced
282 * inadmissible states. We could correct the update now by
283 * performing a limiting step (either convex limiting, or flux
284 * corrected transport)... but *meh*, just perform a backward Euler
285 * step:
286 */
287 step(old_state_vector,
288 t,
289 new_state_vector,
290 tau,
291 id_violation_strategy,
292 reinitialize_gmg,
293 /*crank_nicolson_extrapolation = */ false);
294 }
295 }
296
297
298 template <typename Description, int dim, typename Number>
300 const StateVector &old_state_vector,
301 const Number t,
302 StateVector &new_state_vector,
303 Number tau,
304 const IDViolationStrategy id_violation_strategy,
305 const bool reinitialize_gmg,
306 const bool crank_nicolson_extrapolation) const
307 {
308#ifdef DEBUG_OUTPUT
309 std::cout << "ParabolicSolver<dim, Number>::step()" << std::endl;
310#endif
311
312 constexpr ScalarNumber eps = std::numeric_limits<ScalarNumber>::epsilon();
313
314 const auto &old_U = std::get<0>(old_state_vector);
315 auto &new_U = std::get<0>(new_state_vector);
316
318
319 using VA = VectorizedArray<Number>;
320
321 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
322 const auto &affine_constraints = offline_data_->affine_constraints();
323
324 /* Index ranges for the iteration over the sparsity pattern : */
325
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;
329
330 const auto &sparsity_simd = offline_data_->sparsity_pattern_simd();
331
332 DiagonalMatrix<dim, Number> diagonal_matrix;
333
334#ifdef DEBUG_OUTPUT
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;
338#endif
339
340 /*
341 * A boolean indicating that a restart is required.
342 *
343 * In our current implementation we set this boolean to true if the
344 * backward Euler step produces an internal energy update that
345 * violates the minimum principle, i.e., the minimum of the new
346 * internal energy is smaller than the minimum of the old internal
347 * energy.
348 *
349 * Depending on the chosen "id_violation_strategy" we either signal a
350 * restart by throwing a "Restart" object, or we simply increase the
351 * number of warnings.
352 */
353 std::atomic<bool> restart_needed = false;
354
355 /*
356 * A boolean indicating that we have to correct the high-order Crank
357 * Nicolson update. Note that this is a truly exceptional case
358 * indicating that the high-order update produced an inadmissible
359 * state, *boo*.
360 *
361 * Our current limiting strategy is to simply fall back to perform a
362 * single backward Euler step...
363 */
364 std::atomic<bool> correction_needed = false;
365
366 /*
367 * Step 1:
368 *
369 * Build right hand side for the velocity update.
370 * Also initialize solution vectors for internal energy and velocity
371 * update.
372 */
373 {
374 Scope scope(computing_timer_, "time step [P] 1 - update velocities");
376 LIKWID_MARKER_START("time_step_parabolic_1");
377
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>;
381
382 const auto view = hyperbolic_system_->template view<dim, T>();
383
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);
391
392 write_entry<T>(density_, rho_i, i);
393 /* (5.4a) */
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);
397 }
398 write_entry<T>(internal_energy_, rho_e_i / rho_i, i);
399 }
400 };
401
402 /* Parallel non-vectorized loop: */
403 loop(Number(), n_regular, n_owned);
404 /* Parallel vectorized SIMD loop: */
405 loop(VA(), 0, n_regular);
406
408
409 /*
410 * Set up "strongly enforced" boundary conditions that are not stored
411 * in the AffineConstraints map. In this case we enforce boundary
412 * values by imposing them strongly in the iteration by setting the
413 * initial vector and the right hand side to the right value:
414 */
415
416 const auto &boundary_map = offline_data_->boundary_map();
417
418 for (auto entry : boundary_map) {
419 // [i, normal, normal_mass, boundary_mass, id, position] = entry
420 const auto i = std::get<0>(entry);
421 if (i >= n_owned)
422 continue;
423
424 const auto normal = std::get<1>(entry);
425 const auto id = std::get<4>(entry);
426 const auto position = std::get<5>(entry);
427
428 if (id == Boundary::slip) {
429 /* Remove normal component of velocity: */
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);
435 }
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];
441 }
442
443 } else if (id == Boundary::no_slip) {
444
445 /* Set velocity to zero: */
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.);
449 }
450
451 } else if (id == Boundary::dirichlet) {
452
453 /* Prescribe velocity: */
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;
459
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];
463 }
464
465 internal_energy_.local_element(i) = e_i;
466 }
467 }
468
469 /*
470 * Zero out constrained degrees of freedom due to hanging nodes and
471 * periodic boundary conditions. These boundary conditions are
472 * enforced by modifying the stencil - consequently we have to
473 * remove constrained dofs from the linear system.
474 */
475
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));
481 }
482
483 /* Prepare preconditioner: */
484
485 diagonal_matrix.reinit(
486 lumped_mass_matrix, density_, affine_constraints);
487
488 /*
489 * Update MG matrices all 4 time steps; this is a balance because more
490 * refreshes will render the approximation better, at some additional
491 * cost.
492 */
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());
500
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_);
505
506 for (unsigned int level = level_matrix_free_.min_level();
507 level <= level_matrix_free_.max_level();
508 ++level) {
509 level_velocity_matrices_[level].initialize(
510 *parabolic_system_,
511 *offline_data_,
512 level_matrix_free_[level],
513 level_density_[level],
514 tau,
515 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;
522 } else {
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_;
529 }
530 }
531 mg_smoother_velocity_.initialize(level_velocity_matrices_,
532 smoother_data);
533 }
534
535 LIKWID_MARKER_STOP("time_step_parabolic_1");
536 }
537
538 Number e_min_old;
539
540 {
541 Scope scope(computing_timer_,
542 "time step [X] _ - synchronization barriers");
543
544 /* Compute the global minimum of the internal energy: */
545
546 // .begin() and .end() denote the locally owned index range:
547 e_min_old =
548 *std::min_element(internal_energy_.begin(), internal_energy_.end());
549
550 e_min_old = Utilities::MPI::min(e_min_old,
551 mpi_ensemble_.ensemble_communicator());
552
553 // FIXME: create a meaningful relaxation based on global mesh size min.
554 constexpr Number eps = std::numeric_limits<Number>::epsilon();
555 e_min_old *= (1. - 1000. * eps);
556 }
557
558 /*
559 * Step 1: Solve velocity update:
560 */
561 {
562 Scope scope(computing_timer_, "time step [P] 1 - update velocities");
563
564 LIKWID_MARKER_START("time_step_parabolic_1");
565
566 VelocityMatrix<dim, Number, Number> velocity_operator;
567 velocity_operator.initialize(
568 *parabolic_system_, *offline_data_, matrix_free_, density_, tau);
569
570 const auto tolerance_velocity =
571 (tolerance_linfty_norm_ ? velocity_rhs_.linfty_norm()
572 : velocity_rhs_.l2_norm()) *
573 tolerance_;
574
575 /*
576 * Multigrid might lack robustness for some cases, so in case it takes
577 * too many iterations we better switch to the more robust plain
578 * conjugate gradient method.
579 */
580 try {
581 if (!use_gmg_velocity_)
582 throw SolverControl::NoConvergence(0, 0.);
583
584 using bvt_float = LinearAlgebra::distributed::BlockVector<float>;
585
586 MGCoarseGridApplySmoother<bvt_float> mg_coarse;
587 mg_coarse.initialize(mg_smoother_velocity_);
588
589 mg::Matrix<bvt_float> mg_matrix(level_velocity_matrices_);
590
591 Multigrid<bvt_float> mg(mg_matrix,
592 mg_coarse,
593 mg_transfer_velocity_,
594 mg_smoother_velocity_,
595 mg_smoother_velocity_,
596 level_velocity_matrices_.min_level(),
597 level_velocity_matrices_.max_level());
598
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_);
602
603 SolverControl solver_control(gmg_max_iter_vel_, tolerance_velocity);
604 SolverCG<BlockVector> solver(solver_control);
605 solver.solve(
606 velocity_operator, velocity_, velocity_rhs_, preconditioner);
607
608 /* update exponential moving average */
609 n_iterations_velocity_ =
610 0.9 * n_iterations_velocity_ + 0.1 * solver_control.last_step();
611
612 } catch (SolverControl::NoConvergence &) {
613
614 SolverControl solver_control(1000, tolerance_velocity);
615 SolverCG<BlockVector> solver(solver_control);
616 solver.solve(
617 velocity_operator, velocity_, velocity_rhs_, diagonal_matrix);
618
619 /* update exponential moving average, counting also GMG iterations */
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();
624 }
625
626 LIKWID_MARKER_STOP("time_step_parabolic_1");
627 }
628
629 /*
630 * Step 2: Build internal energy right hand side:
631 */
632 {
633 Scope scope(computing_timer_,
634 "time step [P] 2 - update internal energy");
635
636 LIKWID_MARKER_START("time_step_parabolic_2");
637
638 /* Compute m_i K_i^{n+1/2}: (5.5) */
639 matrix_free_.template cell_loop<ScalarVector, BlockVector>(
640 [this](const auto &data,
641 auto &dst,
642 const auto &src,
643 const auto cell_range) {
644 FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
645 data);
646 FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(data);
647
648 const auto mu = parabolic_system_->mu();
649 const auto lambda = parabolic_system_->lambda();
650
651 for (unsigned int cell = cell_range.first;
652 cell < cell_range.second;
653 ++cell) {
654 velocity.reinit(cell);
655 energy.reinit(cell);
656 velocity.gather_evaluate(src, EvaluationFlags::gradients);
657
658 for (unsigned int q = 0; q < velocity.n_q_points; ++q) {
659 if constexpr (dim == 1) {
660 /* Workaround: no symmetric gradient for 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);
664
665 } else {
666
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);
674 }
675 }
676 energy.integrate_scatter(EvaluationFlags::values, dst);
677 }
678 },
679 internal_energy_rhs_,
680 velocity_,
681 /* zero destination */ true);
682
683 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
684
686
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>;
690
691 const auto view = hyperbolic_system_->template view<dim, T>();
692
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);
699
700 const auto U_i = old_U.template get_tensor<T>(i);
701 const auto V_i = view.momentum(U_i) / rho_i;
702
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);
706 }
707
708 /*
709 * For backward Euler we have to add this algebraic correction
710 * to ensure conservation of total energy.
711 */
712 const auto correction =
713 crank_nicolson_extrapolation
714 ? T(0.)
715 : Number(0.5) * (V_i - V_i_new).norm_square();
716
717 /* rhs_i contains already m_i K_i^{n+1/2} */
718 const auto result = m_i * rho_i * (e_i + correction) + tau * rhs_i;
719 write_entry<T>(internal_energy_rhs_, result, i);
720 }
721 };
722
723 /* Parallel non-vectorized loop: */
724 loop(Number(), n_regular, n_owned);
725 /* Parallel vectorized SIMD loop: */
726 loop(VA(), 0, n_regular);
727
729
730 /*
731 * Set up "strongly enforced" boundary conditions that are not stored
732 * in the AffineConstraints map: We enforce Neumann conditions (i.e.,
733 * insulating boundary conditions) everywhere except for Dirichlet
734 * boundaries where we have to enforce prescribed conditions:
735 */
736
737 const auto &boundary_map = offline_data_->boundary_map();
738
739 for (auto entry : boundary_map) {
740 // [i, normal, normal_mass, boundary_mass, id, position] = entry
741 const auto i = std::get<0>(entry);
742 if (i >= n_owned)
743 continue;
744
745 const auto id = std::get<4>(entry);
746 const auto position = std::get<5>(entry);
747
748 if (id == Boundary::dirichlet) {
749 /* Prescribe internal energy: */
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;
755 }
756 }
757
758 /*
759 * Zero out constrained degrees of freedom due to hanging nodes and
760 * periodic boundary conditions. These boundary conditions are
761 * enforced by modifying the stencil - consequently we have to
762 * remove constrained dofs from the linear system.
763 */
764 affine_constraints.set_zero(internal_energy_);
765 affine_constraints.set_zero(internal_energy_rhs_);
766
767 /*
768 * Update MG matrices all 4 time steps; this is a balance because more
769 * refreshes will render the approximation better, at some additional
770 * cost.
771 */
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());
778
779 level_energy_matrices_.resize(level_matrix_free_.min_level(),
780 level_matrix_free_.max_level());
781
782 for (unsigned int level = level_matrix_free_.min_level();
783 level <= level_matrix_free_.max_level();
784 ++level) {
785 level_energy_matrices_[level].initialize(
786 *offline_data_,
787 level_matrix_free_[level],
788 level_density_[level],
789 tau * parabolic_system_->cv_inverse_kappa(),
790 level);
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;
797 } else {
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_;
804 }
805 }
806 mg_smoother_energy_.initialize(level_energy_matrices_, smoother_data);
807 }
808
809 LIKWID_MARKER_STOP("time_step_parabolic_2");
810 }
811
812 /*
813 * Step 2: Solve internal energy update:
814 */
815 {
816 Scope scope(computing_timer_,
817 "time step [P] 2 - update internal energy");
818
819 LIKWID_MARKER_START("time_step_parabolic_2");
820
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);
825
826 const auto tolerance_internal_energy =
827 (tolerance_linfty_norm_ ? internal_energy_rhs_.linfty_norm()
828 : internal_energy_rhs_.l2_norm()) *
829 tolerance_;
830
831 try {
832 if (!use_gmg_internal_energy_)
833 throw SolverControl::NoConvergence(0, 0.);
834
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_);
839
840 Multigrid<vt_float> mg(mg_matrix,
841 mg_coarse,
842 mg_transfer_energy_,
843 mg_smoother_energy_,
844 mg_smoother_energy_,
845 level_energy_matrices_.min_level(),
846 level_energy_matrices_.max_level());
847
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_);
851
852 SolverControl solver_control(gmg_max_iter_en_,
853 tolerance_internal_energy);
854 SolverCG<ScalarVector> solver(solver_control);
855 solver.solve(energy_operator,
856 internal_energy_,
857 internal_energy_rhs_,
858 preconditioner);
859
860 /* update exponential moving average */
861 n_iterations_internal_energy_ = 0.9 * n_iterations_internal_energy_ +
862 0.1 * solver_control.last_step();
863
864 } catch (SolverControl::NoConvergence &) {
865
866 SolverControl solver_control(1000, tolerance_internal_energy);
867 SolverCG<ScalarVector> solver(solver_control);
868 solver.solve(energy_operator,
869 internal_energy_,
870 internal_energy_rhs_,
871 diagonal_matrix);
872
873 /* update exponential moving average, counting also GMG iterations */
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();
878 }
879
880 LIKWID_MARKER_STOP("time_step_parabolic_2");
881 }
882
883 /*
884 * Step 3: Copy vectors and check for local minimum principle on
885 * internal energy:
886 *
887 * FIXME: Memory access is suboptimal...
888 */
889 {
890 Scope scope(computing_timer_, "time step [P] 3 - write back vectors");
891
893 LIKWID_MARKER_START("time_step_parabolic_3");
894
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>;
898
899 const auto view = hyperbolic_system_->template view<dim, T>();
900
902 for (unsigned int i = left; i < right; i += stride_size) {
903
904 /* Skip constrained degrees of freedom: */
905 const unsigned int row_length = sparsity_simd.row_length(i);
906 if (row_length == 1)
907 continue;
908
909 auto U_i = old_U.template get_tensor<T>(i);
910 const auto rho_i = view.density(U_i);
911
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);
915 }
916
917 auto rho_e_i_new = rho_i * get_entry<T>(internal_energy_, i);
918
919 /*
920 * Check that the backward Euler step itself (which is our "low
921 * order" update) satisfies bounds. If not, signal a restart.
922 */
923
924 if (!(T(0.) == std::max(T(0.), rho_i * e_min_old - rho_e_i_new))) {
925#ifdef DEBUG_OUTPUT
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): "
931 << negative_part(e_i_new - e_min_old) << "\n"
932 << "\t\te_min_new: " << e_i_new << "\n"
933 << std::endl;
934#endif
935 restart_needed = true;
936 }
937
938 if (crank_nicolson_extrapolation) {
939 m_i_new = Number(2.0) * m_i_new - view.momentum(U_i);
940 rho_e_i_new =
941 Number(2.0) * rho_e_i_new - view.internal_energy(U_i);
942
943 /*
944 * If we do perform an extrapolation step for Crank Nicolson
945 * we have to check whether we maintain admissibility
946 */
947
948 if (!(T(0.) ==
949 std::max(T(0.), eps * rho_i * e_min_old - rho_e_i_new))) {
950#ifdef DEBUG_OUTPUT
951 std::cout << std::fixed << std::setprecision(16);
952 const auto e_i_new = rho_e_i_new / rho_i;
953
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;
957#endif
958 correction_needed = true;
959 }
960 }
961
962 const auto E_i_new = rho_e_i_new + 0.5 * m_i_new * m_i_new / rho_i;
963
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;
967
968 new_U.template write_tensor<T>(U_i, i);
969 }
970 };
971
972 /* Parallel non-vectorized loop: */
973 loop(Number(), n_regular, n_owned);
974 /* Parallel vectorized SIMD loop: */
975 loop(VA(), 0, n_regular);
976
977 LIKWID_MARKER_STOP("time_step_parabolic_3");
979
980 new_U.update_ghost_values();
981 }
982
984
985 {
986 Scope scope(computing_timer_,
987 "time step [X] _ - synchronization barriers");
988
989 /*
990 * Synchronize whether we have to restart or correct the time step.
991 * Even though the restart/correction condition itself only affects
992 * the local ensemble we nevertheless need to synchronize the
993 * boolean in case we perform synchronized global time steps.
994 * (Otherwise different ensembles might end up with a different
995 * time step.)
996 */
997
998 restart_needed.store(Utilities::MPI::logical_or(
999 restart_needed.load(),
1000 mpi_ensemble_.synchronization_communicator()));
1001
1002 correction_needed.store(Utilities::MPI::logical_or(
1003 correction_needed.load(),
1004 mpi_ensemble_.synchronization_communicator()));
1005 }
1006
1007 if (correction_needed) {
1008 /* If we can do a restart try that first: */
1009 if (id_violation_strategy == IDViolationStrategy::raise_exception) {
1010 n_restarts_++;
1011 /* Half step size is a good heuristic: */
1012 throw Restart{Number(0.5) * tau};
1013 } else {
1014 n_corrections_++;
1015 throw Correction();
1016 }
1017 }
1018
1019 if (restart_needed) {
1020 switch (id_violation_strategy) {
1022 n_warnings_++;
1023 break;
1025 n_restarts_++;
1026 /* Half step size is a good heuristic: */
1027 throw Restart{Number(0.5) * tau};
1028 }
1029 }
1030 }
1031
1032
1033 template <typename Description, int dim, typename Number>
1035 std::ostream &output) const
1036 {
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 ]")
1042 << std::endl;
1043 }
1044
1045 } // namespace NavierStokes
1046} /* namespace ryujin */
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 &parabolic_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)
Definition: geometry_step.h:25
#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:125
#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
Definition: state_vector.h:53