ryujin 2.1.1 revision 6dc06e5864abd5d99e5d7ab641dbe621936411d9
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 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 NavierStokes module currently only "
135 "supports cG Q1 finite elements."));
136
137 /* Initialize vectors: */
138
139 typename MatrixFree<dim, Number>::AdditionalData additional_data;
140 additional_data.tasks_parallel_scheme =
141 MatrixFree<dim, Number>::AdditionalData::none;
142
143 matrix_free_.reinit(discretization.mapping(),
144 offline_data_->dof_handler(),
145 offline_data_->affine_constraints(),
146 discretization.quadrature_1d(),
147 additional_data);
148
149 const auto &scalar_partitioner =
150 matrix_free_.get_dof_info(0).vector_partitioner;
151
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);
157 }
158
159 internal_energy_.reinit(scalar_partitioner);
160 internal_energy_rhs_.reinit(scalar_partitioner);
161
162 density_.reinit(scalar_partitioner);
163
164 /* Initialize multigrid: */
165
166 if (!use_gmg_velocity_ && !use_gmg_internal_energy_)
167 return;
168
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(),
177 relevant_sets);
178 std::set<types::boundary_id> boundary_ids;
179 boundary_ids.insert(Boundary::dirichlet);
180 boundary_ids.insert(Boundary::no_slip);
181 mg_constrained_dofs_.make_zero_boundary_constraints(
182 offline_data_->dof_handler(), boundary_ids);
183
184 typename MatrixFree<dim, float>::AdditionalData additional_data_level;
185 additional_data_level.tasks_parallel_scheme =
186 MatrixFree<dim, float>::AdditionalData::none;
187
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]);
193 // constraints.add_lines(mg_constrained_dofs_.get_boundary_indices(level));
194 // constraints.merge(mg_constrained_dofs_.get_level_constraints(level));
195 constraints.close();
196 level_matrix_free_[level].reinit(discretization.mapping(),
197 offline_data_->dof_handler(),
198 constraints,
199 discretization.quadrature_1d(),
200 additional_data_level);
201 level_matrix_free_[level].initialize_dof_vector(level_density_[level]);
202 }
203
204 mg_transfer_velocity_.build(offline_data_->dof_handler(),
205 mg_constrained_dofs_,
206 level_matrix_free_);
207 mg_transfer_energy_.build(offline_data_->dof_handler(),
208 level_matrix_free_);
209 }
210
211
212 template <typename Description, int dim, typename Number>
214 const StateVector &old_state_vector,
215 const Number t,
216 StateVector &new_state_vector,
217 Number tau,
218 const IDViolationStrategy id_violation_strategy,
219 const bool reinitialize_gmg) const
220 {
221 /* Backward Euler step to half time step, and extrapolate: */
222
223 step(old_state_vector,
224 t,
225 new_state_vector,
226 tau,
227 id_violation_strategy,
228 reinitialize_gmg,
229 /*crank_nicolson_extrapolation = */ false);
230 }
231
232
233 template <typename Description, int dim, typename Number>
235 const StateVector &old_state_vector,
236 const Number t,
237 StateVector &new_state_vector,
238 Number tau,
239 const IDViolationStrategy id_violation_strategy,
240 const bool reinitialize_gmg) const
241 {
242 try {
243 step(old_state_vector,
244 t,
245 new_state_vector,
246 tau / Number(2.),
247 id_violation_strategy,
248 reinitialize_gmg,
249 /*crank_nicolson_extrapolation = */ true);
250
251 } catch (Correction) {
252
253 /*
254 * Under very rare circumstances we might fail to perform a Crank
255 * Nicolson step because the extrapolation step produced
256 * inadmissible states. We could correct the update now by
257 * performing a limiting step (either convex limiting, or flux
258 * corrected transport)... but *meh*, just perform a backward Euler
259 * step:
260 */
261 step(old_state_vector,
262 t,
263 new_state_vector,
264 tau,
265 id_violation_strategy,
266 reinitialize_gmg,
267 /*crank_nicolson_extrapolation = */ false);
268 }
269 }
270
271
272 template <typename Description, int dim, typename Number>
274 const StateVector &old_state_vector,
275 const Number t,
276 StateVector &new_state_vector,
277 Number tau,
278 const IDViolationStrategy id_violation_strategy,
279 const bool reinitialize_gmg,
280 const bool crank_nicolson_extrapolation) const
281 {
282#ifdef DEBUG_OUTPUT
283 std::cout << "ParabolicSolver<dim, Number>::step()" << std::endl;
284#endif
285
286 constexpr ScalarNumber eps = std::numeric_limits<ScalarNumber>::epsilon();
287
288 const auto &old_U = std::get<0>(old_state_vector);
289 auto &new_U = std::get<0>(new_state_vector);
290
292
293 using VA = VectorizedArray<Number>;
294
295 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
296 const auto &affine_constraints = offline_data_->affine_constraints();
297
298 /* Index ranges for the iteration over the sparsity pattern : */
299
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;
303
304 const auto &sparsity_simd = offline_data_->sparsity_pattern_simd();
305
306 DiagonalMatrix<dim, Number> diagonal_matrix;
307
308#ifdef DEBUG_OUTPUT
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;
312#endif
313
314 /*
315 * A boolean indicating that a restart is required.
316 *
317 * In our current implementation we set this boolean to true if the
318 * backward Euler step produces an internal energy update that
319 * violates the minimum principle, i.e., the minimum of the new
320 * internal energy is smaller than the minimum of the old internal
321 * energy.
322 *
323 * Depending on the chosen "id_violation_strategy" we either signal a
324 * restart by throwing a "Restart" object, or we simply increase the
325 * number of warnings.
326 */
327 std::atomic<bool> restart_needed = false;
328
329 /*
330 * A boolean indicating that we have to correct the high-order Crank
331 * Nicolson update. Note that this is a truly exceptional case
332 * indicating that the high-order update produced an inadmissible
333 * state, *boo*.
334 *
335 * Our current limiting strategy is to simply fall back to perform a
336 * single backward Euler step...
337 */
338 std::atomic<bool> correction_needed = false;
339
340 /*
341 * Step 1:
342 *
343 * Build right hand side for the velocity update.
344 * Also initialize solution vectors for internal energy and velocity
345 * update.
346 */
347 {
348 Scope scope(computing_timer_, "time step [P] 1 - update velocities");
350 LIKWID_MARKER_START("time_step_parabolic_1");
351
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>;
355
356 const auto view = hyperbolic_system_->template view<dim, T>();
357
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);
365
366 write_entry<T>(density_, rho_i, i);
367 /* (5.4a) */
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);
371 }
372 write_entry<T>(internal_energy_, rho_e_i / rho_i, i);
373 }
374 };
375
376 /* Parallel non-vectorized loop: */
377 loop(Number(), n_regular, n_owned);
378 /* Parallel vectorized SIMD loop: */
379 loop(VA(), 0, n_regular);
380
382
383 /*
384 * Set up "strongly enforced" boundary conditions that are not stored
385 * in the AffineConstraints map. In this case we enforce boundary
386 * values by imposing them strongly in the iteration by setting the
387 * initial vector and the right hand side to the right value:
388 */
389
390 const auto &boundary_map = offline_data_->boundary_map();
391
392 for (auto entry : boundary_map) {
393 // [i, normal, normal_mass, boundary_mass, id, position] = entry
394 const auto i = std::get<0>(entry);
395 if (i >= n_owned)
396 continue;
397
398 const auto normal = std::get<1>(entry);
399 const auto id = std::get<4>(entry);
400 const auto position = std::get<5>(entry);
401
402 if (id == Boundary::slip) {
403 /* Remove normal component of velocity: */
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);
409 }
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];
415 }
416
417 } else if (id == Boundary::no_slip) {
418
419 /* Set velocity to zero: */
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.);
423 }
424
425 } else if (id == Boundary::dirichlet) {
426
427 /* Prescribe velocity: */
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;
433
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];
437 }
438
439 internal_energy_.local_element(i) = e_i;
440 }
441 }
442
443 /*
444 * Zero out constrained degrees of freedom due to hanging nodes and
445 * periodic boundary conditions. These boundary conditions are
446 * enforced by modifying the stencil - consequently we have to
447 * remove constrained dofs from the linear system.
448 */
449
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));
455 }
456
457 /* Prepare preconditioner: */
458
459 diagonal_matrix.reinit(
460 lumped_mass_matrix, density_, affine_constraints);
461
462 /*
463 * Update MG matrices all 4 time steps; this is a balance because more
464 * refreshes will render the approximation better, at some additional
465 * cost.
466 */
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());
474
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_);
479
480 for (unsigned int level = level_matrix_free_.min_level();
481 level <= level_matrix_free_.max_level();
482 ++level) {
483 level_velocity_matrices_[level].initialize(
484 *parabolic_system_,
485 *offline_data_,
486 level_matrix_free_[level],
487 level_density_[level],
488 tau,
489 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;
496 } else {
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_;
503 }
504 }
505 mg_smoother_velocity_.initialize(level_velocity_matrices_,
506 smoother_data);
507 }
508
509 LIKWID_MARKER_STOP("time_step_parabolic_1");
510 }
511
512 Number e_min_old;
513
514 {
515 Scope scope(computing_timer_,
516 "time step [P] _ - synchronization barriers");
517
518 /* Compute the global minimum of the internal energy: */
519
520 // .begin() and .end() denote the locally owned index range:
521 e_min_old =
522 *std::min_element(internal_energy_.begin(), internal_energy_.end());
523
524 e_min_old = Utilities::MPI::min(e_min_old,
525 mpi_ensemble_.ensemble_communicator());
526
527 // FIXME: create a meaningful relaxation based on global mesh size min.
528 constexpr Number eps = std::numeric_limits<Number>::epsilon();
529 e_min_old *= (1. - 1000. * eps);
530 }
531
532 /*
533 * Step 1: Solve velocity update:
534 */
535 {
536 Scope scope(computing_timer_, "time step [P] 1 - update velocities");
537
538 LIKWID_MARKER_START("time_step_parabolic_1");
539
540 VelocityMatrix<dim, Number, Number> velocity_operator;
541 velocity_operator.initialize(
542 *parabolic_system_, *offline_data_, matrix_free_, density_, tau);
543
544 const auto tolerance_velocity =
545 (tolerance_linfty_norm_ ? velocity_rhs_.linfty_norm()
546 : velocity_rhs_.l2_norm()) *
547 tolerance_;
548
549 /*
550 * Multigrid might lack robustness for some cases, so in case it takes
551 * too many iterations we better switch to the more robust plain
552 * conjugate gradient method.
553 */
554 try {
555 if (!use_gmg_velocity_)
556 throw SolverControl::NoConvergence(0, 0.);
557
558 using bvt_float = LinearAlgebra::distributed::BlockVector<float>;
559
560 MGCoarseGridApplySmoother<bvt_float> mg_coarse;
561 mg_coarse.initialize(mg_smoother_velocity_);
562
563 mg::Matrix<bvt_float> mg_matrix(level_velocity_matrices_);
564
565 Multigrid<bvt_float> mg(mg_matrix,
566 mg_coarse,
567 mg_transfer_velocity_,
568 mg_smoother_velocity_,
569 mg_smoother_velocity_,
570 level_velocity_matrices_.min_level(),
571 level_velocity_matrices_.max_level());
572
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_);
576
577 SolverControl solver_control(gmg_max_iter_vel_, tolerance_velocity);
578 SolverCG<BlockVector> solver(solver_control);
579 solver.solve(
580 velocity_operator, velocity_, velocity_rhs_, preconditioner);
581
582 /* update exponential moving average */
583 n_iterations_velocity_ =
584 0.9 * n_iterations_velocity_ + 0.1 * solver_control.last_step();
585
586 } catch (SolverControl::NoConvergence &) {
587
588 SolverControl solver_control(1000, tolerance_velocity);
589 SolverCG<BlockVector> solver(solver_control);
590 solver.solve(
591 velocity_operator, velocity_, velocity_rhs_, diagonal_matrix);
592
593 /* update exponential moving average, counting also GMG iterations */
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();
598 }
599
600 LIKWID_MARKER_STOP("time_step_parabolic_1");
601 }
602
603 /*
604 * Step 2: Build internal energy right hand side:
605 */
606 {
607 Scope scope(computing_timer_,
608 "time step [P] 2 - update internal energy");
609
610 LIKWID_MARKER_START("time_step_parabolic_2");
611
612 /* Compute m_i K_i^{n+1/2}: (5.5) */
613 matrix_free_.template cell_loop<ScalarVector, BlockVector>(
614 [this](const auto &data,
615 auto &dst,
616 const auto &src,
617 const auto cell_range) {
618 FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
619 data);
620 FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(data);
621
622 const auto mu = parabolic_system_->mu();
623 const auto lambda = parabolic_system_->lambda();
624
625 for (unsigned int cell = cell_range.first;
626 cell < cell_range.second;
627 ++cell) {
628 velocity.reinit(cell);
629 energy.reinit(cell);
630 velocity.gather_evaluate(src, EvaluationFlags::gradients);
631
632 for (unsigned int q = 0; q < velocity.n_q_points; ++q) {
633 if constexpr (dim == 1) {
634 /* Workaround: no symmetric gradient for 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);
638
639 } else {
640
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);
648 }
649 }
650 energy.integrate_scatter(EvaluationFlags::values, dst);
651 }
652 },
653 internal_energy_rhs_,
654 velocity_,
655 /* zero destination */ true);
656
657 const auto &lumped_mass_matrix = offline_data_->lumped_mass_matrix();
658
660
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>;
664
665 const auto view = hyperbolic_system_->template view<dim, T>();
666
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);
673
674 const auto U_i = old_U.template get_tensor<T>(i);
675 const auto V_i = view.momentum(U_i) / rho_i;
676
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);
680 }
681
682 /*
683 * For backward Euler we have to add this algebraic correction
684 * to ensure conservation of total energy.
685 */
686 const auto correction =
687 crank_nicolson_extrapolation
688 ? T(0.)
689 : Number(0.5) * (V_i - V_i_new).norm_square();
690
691 /* rhs_i contains already m_i K_i^{n+1/2} */
692 const auto result = m_i * rho_i * (e_i + correction) + tau * rhs_i;
693 write_entry<T>(internal_energy_rhs_, result, i);
694 }
695 };
696
697 /* Parallel non-vectorized loop: */
698 loop(Number(), n_regular, n_owned);
699 /* Parallel vectorized SIMD loop: */
700 loop(VA(), 0, n_regular);
701
703
704 /*
705 * Set up "strongly enforced" boundary conditions that are not stored
706 * in the AffineConstraints map: We enforce Neumann conditions (i.e.,
707 * insulating boundary conditions) everywhere except for Dirichlet
708 * boundaries where we have to enforce prescribed conditions:
709 */
710
711 const auto &boundary_map = offline_data_->boundary_map();
712
713 for (auto entry : boundary_map) {
714 // [i, normal, normal_mass, boundary_mass, id, position] = entry
715 const auto i = std::get<0>(entry);
716 if (i >= n_owned)
717 continue;
718
719 const auto id = std::get<4>(entry);
720 const auto position = std::get<5>(entry);
721
722 if (id == Boundary::dirichlet) {
723 /* Prescribe internal energy: */
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;
729 }
730 }
731
732 /*
733 * Zero out constrained degrees of freedom due to hanging nodes and
734 * periodic boundary conditions. These boundary conditions are
735 * enforced by modifying the stencil - consequently we have to
736 * remove constrained dofs from the linear system.
737 */
738 affine_constraints.set_zero(internal_energy_);
739 affine_constraints.set_zero(internal_energy_rhs_);
740
741 /*
742 * Update MG matrices all 4 time steps; this is a balance because more
743 * refreshes will render the approximation better, at some additional
744 * cost.
745 */
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());
752
753 level_energy_matrices_.resize(level_matrix_free_.min_level(),
754 level_matrix_free_.max_level());
755
756 for (unsigned int level = level_matrix_free_.min_level();
757 level <= level_matrix_free_.max_level();
758 ++level) {
759 level_energy_matrices_[level].initialize(
760 *offline_data_,
761 level_matrix_free_[level],
762 level_density_[level],
763 tau * parabolic_system_->cv_inverse_kappa(),
764 level);
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;
771 } else {
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_;
778 }
779 }
780 mg_smoother_energy_.initialize(level_energy_matrices_, smoother_data);
781 }
782
783 LIKWID_MARKER_STOP("time_step_parabolic_2");
784 }
785
786 /*
787 * Step 2: Solve internal energy update:
788 */
789 {
790 Scope scope(computing_timer_,
791 "time step [P] 2 - update internal energy");
792
793 LIKWID_MARKER_START("time_step_parabolic_2");
794
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);
799
800 const auto tolerance_internal_energy =
801 (tolerance_linfty_norm_ ? internal_energy_rhs_.linfty_norm()
802 : internal_energy_rhs_.l2_norm()) *
803 tolerance_;
804
805 try {
806 if (!use_gmg_internal_energy_)
807 throw SolverControl::NoConvergence(0, 0.);
808
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_);
813
814 Multigrid<vt_float> mg(mg_matrix,
815 mg_coarse,
816 mg_transfer_energy_,
817 mg_smoother_energy_,
818 mg_smoother_energy_,
819 level_energy_matrices_.min_level(),
820 level_energy_matrices_.max_level());
821
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_);
825
826 SolverControl solver_control(gmg_max_iter_en_,
827 tolerance_internal_energy);
828 SolverCG<ScalarVector> solver(solver_control);
829 solver.solve(energy_operator,
830 internal_energy_,
831 internal_energy_rhs_,
832 preconditioner);
833
834 /* update exponential moving average */
835 n_iterations_internal_energy_ = 0.9 * n_iterations_internal_energy_ +
836 0.1 * solver_control.last_step();
837
838 } catch (SolverControl::NoConvergence &) {
839
840 SolverControl solver_control(1000, tolerance_internal_energy);
841 SolverCG<ScalarVector> solver(solver_control);
842 solver.solve(energy_operator,
843 internal_energy_,
844 internal_energy_rhs_,
845 diagonal_matrix);
846
847 /* update exponential moving average, counting also GMG iterations */
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();
852 }
853
854 LIKWID_MARKER_STOP("time_step_parabolic_2");
855 }
856
857 /*
858 * Step 3: Copy vectors and check for local minimum principle on
859 * internal energy:
860 *
861 * FIXME: Memory access is suboptimal...
862 */
863 {
864 Scope scope(computing_timer_, "time step [P] 3 - write back vectors");
865
867 LIKWID_MARKER_START("time_step_parabolic_3");
868
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>;
872
873 const auto view = hyperbolic_system_->template view<dim, T>();
874
876 for (unsigned int i = left; i < right; i += stride_size) {
877
878 /* Skip constrained degrees of freedom: */
879 const unsigned int row_length = sparsity_simd.row_length(i);
880 if (row_length == 1)
881 continue;
882
883 auto U_i = old_U.template get_tensor<T>(i);
884 const auto rho_i = view.density(U_i);
885
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);
889 }
890
891 auto rho_e_i_new = rho_i * get_entry<T>(internal_energy_, i);
892
893 /*
894 * Check that the backward Euler step itself (which is our "low
895 * order" update) satisfies bounds. If not, signal a restart.
896 */
897
898 if (!(T(0.) == std::max(T(0.), rho_i * e_min_old - rho_e_i_new))) {
899#ifdef DEBUG_OUTPUT
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): "
905 << negative_part(e_i_new - e_min_old) << "\n"
906 << "\t\te_min_new: " << e_i_new << "\n"
907 << std::endl;
908#endif
909 restart_needed = true;
910 }
911
912 if (crank_nicolson_extrapolation) {
913 m_i_new = Number(2.0) * m_i_new - view.momentum(U_i);
914 rho_e_i_new =
915 Number(2.0) * rho_e_i_new - view.internal_energy(U_i);
916
917 /*
918 * If we do perform an extrapolation step for Crank Nicolson
919 * we have to check whether we maintain admissibility
920 */
921
922 if (!(T(0.) ==
923 std::max(T(0.), eps * rho_i * e_min_old - rho_e_i_new))) {
924#ifdef DEBUG_OUTPUT
925 std::cout << std::fixed << std::setprecision(16);
926 const auto e_i_new = rho_e_i_new / rho_i;
927
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;
931#endif
932 correction_needed = true;
933 }
934 }
935
936 const auto E_i_new = rho_e_i_new + 0.5 * m_i_new * m_i_new / rho_i;
937
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;
941
942 new_U.template write_tensor<T>(U_i, i);
943 }
944 };
945
946 /* Parallel non-vectorized loop: */
947 loop(Number(), n_regular, n_owned);
948 /* Parallel vectorized SIMD loop: */
949 loop(VA(), 0, n_regular);
950
951 LIKWID_MARKER_STOP("time_step_parabolic_3");
953
954 new_U.update_ghost_values();
955 }
956
958
959 {
960 Scope scope(computing_timer_,
961 "time step [H] _ - synchronization barriers");
962
963 /*
964 * Synchronize whether we have to restart or correct the time step.
965 * Even though the restart/correction condition itself only affects
966 * the local ensemble we nevertheless need to synchronize the
967 * boolean in case we perform synchronized global time steps.
968 * (Otherwise different ensembles might end up with a different
969 * time step.)
970 */
971
972 restart_needed.store(Utilities::MPI::logical_or(
973 restart_needed.load(),
974 mpi_ensemble_.synchronization_communicator()));
975
976 correction_needed.store(Utilities::MPI::logical_or(
977 correction_needed.load(),
978 mpi_ensemble_.synchronization_communicator()));
979 }
980
981 if (correction_needed) {
982 /* If we can do a restart try that first: */
983 if (id_violation_strategy == IDViolationStrategy::raise_exception) {
984 n_restarts_++;
985 throw Restart();
986 } else {
987 n_corrections_++;
988 throw Correction();
989 }
990 }
991
992 if (restart_needed) {
993 switch (id_violation_strategy) {
995 n_warnings_++;
996 break;
998 n_restarts_++;
999 throw Restart();
1000 }
1001 }
1002 }
1003
1004
1005 template <typename Description, int dim, typename Number>
1007 std::ostream &output) const
1008 {
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 ]")
1014 << std::endl;
1015 }
1016
1017 } // namespace NavierStokes
1018} /* 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
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:23
#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
std::tuple< MultiComponentVector< Number, problem_dim >, MultiComponentVector< Number, prec_dim >, BlockVector< Number > > StateVector
Definition: state_vector.h:51