ryujin 2.1.1 revision 46bf70e400e423a8ffffe8300887eeb35b8dfb2c
parabolic_solver_gmg_operators.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 <compile_time_options.h>
9
10#include <observer_pointer.h>
11#include <offline_data.h>
12#include <openmp.h>
13#include <simd.h>
14
15#include "../euler/hyperbolic_system.h"
16#include "parabolic_system.h"
17
18#include <deal.II/base/vectorization.h>
19#include <deal.II/lac/diagonal_matrix.h>
20#include <deal.II/matrix_free/fe_evaluation.h>
21#include <deal.II/multigrid/mg_base.h>
22#include <deal.II/multigrid/mg_transfer_matrix_free.h>
23
24/*
25 * FIXME: generalize and make these operators equation independent and
26 * refactor into ../parabolic_module_gmg_operators.h
27 */
28
29namespace ryujin
30{
31 namespace NavierStokes
32 {
34
43 template <int dim, typename Number>
45 {
46 public:
50 using vector_type = dealii::LinearAlgebra::distributed::Vector<Number>;
51
56 dealii::LinearAlgebra::distributed::BlockVector<Number>;
57
61 DiagonalMatrix() = default;
62
67 void reinit(const vector_type &lumped_mass_matrix,
68 const vector_type &density,
69 const dealii::AffineConstraints<Number> &affine_constraints)
70 {
71 diagonal.reinit(density, true);
72
73 DEAL_II_OPENMP_SIMD_PRAGMA
74 for (unsigned int i = 0;
75 i < density.get_partitioner()->locally_owned_size();
76 ++i) {
77 diagonal.local_element(i) =
78 Number(1.0) /
79 (density.local_element(i) * lumped_mass_matrix.local_element(i));
80 }
81
82 /*
83 * Fix up diagonal entries for constrained degrees of freedom due to
84 * periodic boundary conditions.
85 */
86 affine_constraints.set_zero(diagonal);
87 }
88
93 {
94 return diagonal;
95 }
96
101 {
102 return diagonal_block;
103 }
104
108 void vmult(vector_type &dst, const vector_type &src) const
109 {
110 AssertDimension(diagonal_block.size(), 0);
111 DEAL_II_OPENMP_SIMD_PRAGMA
112 for (unsigned int i = 0;
113 i < diagonal.get_partitioner()->locally_owned_size();
114 ++i)
115 dst.local_element(i) =
116 diagonal.local_element(i) * src.local_element(i);
117 }
118
122 void vmult(block_vector_type &dst, const block_vector_type &src) const
123 {
124 AssertDimension(dim, dst.n_blocks());
125 AssertDimension(dim, src.n_blocks());
126 if (diagonal_block.size() == 0) {
127 DEAL_II_OPENMP_SIMD_PRAGMA
128 for (unsigned int i = 0;
129 i < diagonal.get_partitioner()->locally_owned_size();
130 ++i)
131 for (unsigned int d = 0; d < dim; ++d)
132 dst.block(d).local_element(i) =
133 diagonal.local_element(i) * src.block(d).local_element(i);
134 } else
135 for (unsigned int d = 0; d < dim; ++d) {
136 DEAL_II_OPENMP_SIMD_PRAGMA
137 for (unsigned int i = 0;
138 i < src.block(d).get_partitioner()->locally_owned_size();
139 ++i)
140 dst.block(d).local_element(i) =
141 diagonal_block.block(d).local_element(i) *
142 src.block(d).local_element(i);
143 }
144 }
145
146 private:
147 vector_type diagonal;
148 block_vector_type diagonal_block;
149 };
150
151
158 template <int dim, typename Number, typename Number2>
160 {
161 public:
162 // FIXME: refactor
163 static constexpr unsigned int order_fe = 1;
164 static constexpr unsigned int order_quad = 2;
165
166 using vector_type = dealii::LinearAlgebra::distributed::Vector<Number>;
168 dealii::LinearAlgebra::distributed::BlockVector<Number>;
169
170 VelocityMatrix() = default;
171
173 const ParabolicSystem &parabolic_system,
174 const OfflineData<dim, Number2> &offline_data,
175 const dealii::MatrixFree<dim, Number> &matrix_free,
176 const dealii::LinearAlgebra::distributed::Vector<Number> &density,
177 const Number theta_x_tau,
178 const unsigned int level = dealii::numbers::invalid_unsigned_int)
179 {
180 parabolic_system_ = &parabolic_system;
181 offline_data_ = &offline_data;
182 matrix_free_ = &matrix_free;
183 density_ = &density;
184 theta_x_tau_ = theta_x_tau;
185 level_ = level;
186 }
187
188 void Tvmult(block_vector_type &dst, const block_vector_type &src) const
189 {
190 vmult(dst, src);
191 }
192
193 void vmult(block_vector_type &dst, const block_vector_type &src) const
194 {
195 /* Apply action of m_i rho_i V_i: */
196
197 using VA = dealii::VectorizedArray<Number>;
198 constexpr auto simd_length = VA::size();
199
200 const vector_type *lumped_mass_matrix = nullptr;
201 if constexpr (std::is_same<Number, Number2>::value) {
202 if constexpr (std::is_same<Number, float>::value) {
203 if (level_ == dealii::numbers::invalid_unsigned_int)
204 lumped_mass_matrix = &offline_data_->lumped_mass_matrix();
205 else
206 lumped_mass_matrix =
207 &offline_data_->level_lumped_mass_matrix()[level_];
208 } else {
209 Assert(level_ == dealii::numbers::invalid_unsigned_int,
210 dealii::ExcInternalError());
211 lumped_mass_matrix = &offline_data_->lumped_mass_matrix();
212 }
213 } else
214 lumped_mass_matrix =
215 &offline_data_->level_lumped_mass_matrix()[level_];
216
217 const unsigned int n_owned =
218 lumped_mass_matrix->get_partitioner()->locally_owned_size();
219 const unsigned int size_regular = n_owned / simd_length * simd_length;
220
222
224 for (unsigned int i = 0; i < size_regular; i += simd_length) {
225 const auto m_i = get_entry<VA>(*lumped_mass_matrix, i);
226 const auto rho_i = get_entry<VA>(*density_, i);
227 for (unsigned int d = 0; d < dim; ++d) {
228 const auto temp = get_entry<VA>(src.block(d), i);
229 write_entry<VA>(dst.block(d), m_i * rho_i * temp, i);
230 }
231 }
232
234
235 for (unsigned int i = size_regular; i < n_owned; ++i) {
236 const auto m_i = lumped_mass_matrix->local_element(i);
237 const auto rho_i = density_->local_element(i);
238
239 for (unsigned int d = 0; d < dim; ++d) {
240 const auto temp = src.block(d).local_element(i);
241 dst.block(d).local_element(i) = m_i * rho_i * temp;
242 }
243 }
244
245 /* Apply action of stress tensor: + theta * \sum_j B_ij V_j: */
246
247 const auto integrator = [this](const auto &data,
248 auto &dst,
249 const auto &src,
250 const auto range) {
251 dealii::FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
252 data);
253
254 for (unsigned int cell = range.first; cell < range.second; ++cell) {
255 velocity.reinit(cell);
256 velocity.read_dof_values(src);
257 apply_local_operator(velocity);
258 velocity.distribute_local_to_global(dst);
259 }
260 };
261
262 matrix_free_->template cell_loop<block_vector_type, block_vector_type>(
263 integrator, dst, src, /* zero destination */ false);
264
265 /* (5.4a) Fix up constrained degrees of freedom: */
266
267 const auto &boundary_map =
268 level_ == dealii::numbers::invalid_unsigned_int
269 ? offline_data_->boundary_map()
270 : offline_data_->level_boundary_map()[level_];
271
272 for (auto entry : boundary_map) {
273 // [i, normal, normal_mass, boundary_mass, id, position] = entry
274 const auto i = std::get<0>(entry);
275 if (i >= n_owned)
276 continue;
277
278 const dealii::Tensor<1, dim, Number> normal = std::get<1>(entry);
279 const auto id = std::get<4>(entry);
280
281 if (id == Boundary::slip) {
282 dealii::Tensor<1, dim, Number> V_i;
283 for (unsigned int d = 0; d < dim; ++d)
284 V_i[d] = dst.block(d).local_element(i);
285
286 /* replace normal component by source */
287 V_i -= 1. * (V_i * normal) * normal;
288 for (unsigned int d = 0; d < dim; ++d) {
289 const auto src_d = src.block(d).local_element(i);
290 V_i += 1. * (src_d * normal[d]) * normal;
291 }
292
293 for (unsigned int d = 0; d < dim; ++d)
294 dst.block(d).local_element(i) = V_i[d];
295
296 } else if (id == Boundary::no_slip || id == Boundary::dirichlet) {
297
298 /* set dst to src vector: */
299 for (unsigned int d = 0; d < dim; ++d)
300 dst.block(d).local_element(i) = src.block(d).local_element(i);
301 }
302 }
303 }
304
306 std::shared_ptr<DiagonalMatrix<dim, Number>> &matrix) const
307 {
308 Assert(level_ != dealii::numbers::invalid_unsigned_int,
309 dealii::ExcNotImplemented());
310 matrix = std::make_shared<DiagonalMatrix<dim, Number>>();
311 block_vector_type &vector = matrix->get_block_vector();
312 vector.reinit(dim);
313 for (unsigned int d = 0; d < dim; ++d)
314 matrix_free_->initialize_dof_vector(vector.block(d));
315 vector.collect_sizes();
316
317 const auto &lumped_mass_matrix =
318 offline_data_->level_lumped_mass_matrix()[level_];
319
320 unsigned int dummy = 0;
321 matrix_free_->template cell_loop<block_vector_type, unsigned int>(
322 [this](
323 const auto &data, auto &dst, const auto &, const auto range) {
324 dealii::FEEvaluation<dim, order_fe, order_quad, dim, Number>
325 velocity(data);
326 dealii::FEEvaluation<dim, order_fe, order_quad, dim, Number>
327 writer(data);
328
329 for (unsigned int cell = range.first; cell < range.second;
330 ++cell) {
331 velocity.reinit(cell);
332 writer.reinit(cell);
333 for (unsigned int i = 0; i < velocity.dofs_per_cell; ++i) {
334 for (unsigned int j = 0; j < velocity.dofs_per_cell; ++j)
335 velocity.begin_dof_values()[j] =
336 dealii::VectorizedArray<Number>();
337 velocity.begin_dof_values()[i] =
338 dealii::make_vectorized_array<Number>(1.);
339 apply_local_operator(velocity);
340 writer.begin_dof_values()[i] = velocity.begin_dof_values()[i];
341 }
342 writer.distribute_local_to_global(dst);
343 }
344 },
345 vector,
346 dummy,
347 /* zero destination */ true);
348
349 const unsigned int n_owned =
350 lumped_mass_matrix.get_partitioner()->locally_owned_size();
351
353
355 for (unsigned int i = 0; i < n_owned; ++i) {
356 const auto m_i = lumped_mass_matrix.local_element(i);
357 const auto rho_i = density_->local_element(i);
358 for (unsigned int d = 0; d < dim; ++d)
359 vector.block(d).local_element(i) =
360 1. / (m_i * rho_i + vector.block(d).local_element(i));
361 }
362
364
365 const auto &boundary_map = offline_data_->level_boundary_map()[level_];
366
367 for (auto entry : boundary_map) {
368 // [i, normal, normal_mass, boundary_mass, id, position] = entry
369 const auto i = std::get<0>(entry);
370 if (i >= n_owned)
371 continue;
372
373 const dealii::Tensor<1, dim, Number> normal = std::get<1>(entry);
374 const auto id = std::get<4>(entry);
375
376 if (id == Boundary::slip) {
377 dealii::Tensor<1, dim, Number> V_i;
378 for (unsigned int d = 0; d < dim; ++d)
379 V_i[d] = vector.block(d).local_element(i);
380
381 /* replace normal component by 1. */
382 V_i -= 1. * (V_i * normal) * normal;
383 for (unsigned int d = 0; d < dim; ++d) {
384 V_i += 1. * (1. * normal[d]) * normal;
385 }
386
387 for (unsigned int d = 0; d < dim; ++d)
388 vector.block(d).local_element(i) = V_i[d];
389
390 } else if (id == Boundary::no_slip || id == Boundary::dirichlet) {
391
392 /* set dst to src vector: */
393 for (unsigned int d = 0; d < dim; ++d)
394 vector.block(d).local_element(i) = 1.;
395 }
396 }
397 }
398
399 private:
400 const ParabolicSystem *parabolic_system_;
401 const OfflineData<dim, Number2> *offline_data_;
402 const dealii::MatrixFree<dim, Number> *matrix_free_;
403 const vector_type *density_;
404 Number theta_x_tau_;
405 unsigned int level_;
406
407 template <typename Evaluator>
408 void apply_local_operator(Evaluator &velocity) const
409 {
410 const Number mu = parabolic_system_->mu();
411 const Number lambda = parabolic_system_->lambda();
412
413 velocity.evaluate(dealii::EvaluationFlags::gradients);
414
415 for (unsigned int q = 0; q < velocity.n_q_points; ++q) {
416 if constexpr (dim == 1) {
417 /* Workaround: no symmetric gradient for dim == 1: */
418 const auto gradient = velocity.get_gradient(q);
419 auto S = (4. / 3. * mu + lambda) * gradient;
420 velocity.submit_gradient(theta_x_tau_ * S, q);
421
422 } else {
423
424 const auto symmetric_gradient = velocity.get_symmetric_gradient(q);
425 const auto divergence = trace(symmetric_gradient);
426 // S = (2 mu nabla^S(v) + (lambda - 2/3*mu) div(v) Id) : nabla phi
427 auto S = 2. * mu * symmetric_gradient;
428 for (unsigned int d = 0; d < dim; ++d)
429 S[d][d] += (lambda - 2. / 3. * mu) * divergence;
430 velocity.submit_symmetric_gradient(theta_x_tau_ * S, q);
431 }
432 }
433
434 velocity.integrate(dealii::EvaluationFlags::gradients);
435 }
436 };
437
438
442 template <int dim, typename Number>
444 : public dealii::MGTransferBase<
445 dealii::LinearAlgebra::distributed::BlockVector<Number>>
446 {
447 public:
448 using scalar_type = dealii::LinearAlgebra::distributed::Vector<Number>;
450 dealii::LinearAlgebra::distributed::BlockVector<Number>;
451
453
454 void build(const dealii::DoFHandler<dim> &dof_handler,
455 const dealii::MGConstrainedDoFs &mg_constrained_dofs,
456 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
457 &matrix_free)
458 {
459 transfer_.initialize_constraints(mg_constrained_dofs);
460 transfer_.build(dof_handler);
461 level_matrix_free_ = &matrix_free;
462 scalar_vector.resize(matrix_free.min_level(), matrix_free.max_level());
463 for (unsigned int level = matrix_free.min_level();
464 level < matrix_free.max_level();
465 ++level)
466 matrix_free[level].initialize_dof_vector(scalar_vector[level]);
467 }
468
469 void prolongate(const unsigned int to_level,
470 vector_type &dst,
471 const vector_type &src) const override
472 {
473 for (unsigned int block = 0; block < src.n_blocks(); ++block)
474 transfer_.prolongate(to_level, dst.block(block), src.block(block));
475 }
476
477 void restrict_and_add(const unsigned int to_level,
478 vector_type &dst,
479 const vector_type &src) const override
480 {
481 for (unsigned int block = 0; block < src.n_blocks(); ++block)
482 transfer_.restrict_and_add(
483 to_level, dst.block(block), src.block(block));
484 }
485
486 template <typename Number2>
488 const dealii::DoFHandler<dim> &dof_handler,
489 dealii::MGLevelObject<scalar_type> &dst,
490 const dealii::LinearAlgebra::distributed::Vector<Number2> &src) const
491 {
492 if (dst[dst.min_level()].size() == 0)
493 for (unsigned int l = dst.min_level(); l <= dst.max_level(); ++l)
494 (*level_matrix_free_)[l].initialize_dof_vector(dst[l]);
495 transfer_.interpolate_to_mg(dof_handler, dst, src);
496 }
497
498 template <typename Number2>
499 void
500 copy_to_mg(const dealii::DoFHandler<dim> &dof_handler,
501 dealii::MGLevelObject<vector_type> &dst,
502 const dealii::LinearAlgebra::distributed::BlockVector<Number2>
503 &src) const
504 {
505 if (dst[dst.min_level()].size() == 0)
506 for (unsigned int l = dst.min_level(); l <= dst.max_level(); ++l) {
507 dst[l].reinit(src.n_blocks());
508 for (unsigned int block = 0; block < src.n_blocks(); ++block)
509 (*level_matrix_free_)[l].initialize_dof_vector(
510 dst[l].block(block));
511 dst[l].collect_sizes();
512 }
513
514 for (unsigned int block = 0; block < src.n_blocks(); ++block) {
515 transfer_.copy_to_mg(dof_handler, scalar_vector, src.block(block));
516 for (unsigned int level = dst.min_level(); level <= dst.max_level();
517 ++level)
518 dst[level].block(block).copy_locally_owned_data_from(
519 scalar_vector[level]);
520 }
521 }
522
523 template <typename Number2>
525 const dealii::DoFHandler<dim> &dof_handler,
526 dealii::LinearAlgebra::distributed::BlockVector<Number2> &dst,
527 const dealii::MGLevelObject<vector_type> &src) const
528 {
529 for (unsigned int block = 0; block < dst.n_blocks(); ++block) {
530 for (unsigned int level = src.min_level(); level <= src.max_level();
531 ++level)
532 scalar_vector[level].copy_locally_owned_data_from(
533 src[level].block(block));
534 transfer_.copy_from_mg(dof_handler, dst.block(block), scalar_vector);
535 }
536 }
537
538 private:
539 dealii::MGTransferMatrixFree<dim, Number> transfer_;
540 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
541 *level_matrix_free_;
542 mutable dealii::MGLevelObject<scalar_type> scalar_vector;
543 };
544
545
549 template <int dim, typename Number, typename Number2>
551 {
552 public:
553 // FIXME: refactor
554 static constexpr unsigned int order_fe = 1;
555 static constexpr unsigned int order_quad = 2;
556
557 using vector_type = dealii::LinearAlgebra::distributed::Vector<Number>;
558
559 EnergyMatrix() = default;
560
562 const OfflineData<dim, Number2> &offline_data,
563 const dealii::MatrixFree<dim, Number> &matrix_free,
564 const dealii::LinearAlgebra::distributed::Vector<Number> &density,
565 const Number time_factor,
566 const unsigned int level = dealii::numbers::invalid_unsigned_int)
567 {
568 offline_data_ = &offline_data;
569 matrix_free_ = &matrix_free;
570 density_ = &density;
571 factor_ = time_factor;
572 level_ = level;
573 }
574
575 void Tvmult(vector_type &dst, const vector_type &src) const
576 {
577 vmult(dst, src);
578 }
579
580 dealii::types::global_dof_index m() const
581 {
582 return density_->size();
583 }
584
585 Number el(const unsigned int, const unsigned int) const
586 {
587 Assert(false, dealii::ExcNotImplemented());
588 return Number();
589 }
590
591 void vmult(vector_type &dst, const vector_type &src) const
592 {
593 /* Apply action of m_i rho_i V_i: */
594
595 using VA = dealii::VectorizedArray<Number>;
596 constexpr auto simd_length = VA::size();
597
598 const vector_type *lumped_mass_matrix = nullptr;
599 if constexpr (std::is_same<Number, Number2>::value) {
600 if constexpr (std::is_same<Number, float>::value) {
601 if (level_ == dealii::numbers::invalid_unsigned_int)
602 lumped_mass_matrix = &offline_data_->lumped_mass_matrix();
603 else
604 lumped_mass_matrix =
605 &offline_data_->level_lumped_mass_matrix()[level_];
606 } else {
607 Assert(level_ == dealii::numbers::invalid_unsigned_int,
608 dealii::ExcInternalError());
609 lumped_mass_matrix = &offline_data_->lumped_mass_matrix();
610 }
611 } else
612 lumped_mass_matrix =
613 &offline_data_->level_lumped_mass_matrix()[level_];
614
615 const unsigned int n_owned =
616 lumped_mass_matrix->get_partitioner()->locally_owned_size();
617 const unsigned int size_regular = n_owned / simd_length * simd_length;
618
620
622 for (unsigned int i = 0; i < size_regular; i += simd_length) {
623 const auto m_i = get_entry<VA>(*lumped_mass_matrix, i);
624 const auto rho_i = get_entry<VA>(*density_, i);
625 const auto e_i = get_entry<VA>(src, i);
626 write_entry<VA>(dst, m_i * rho_i * e_i, i);
627 }
628
630
631 for (unsigned int i = size_regular; i < n_owned; ++i) {
632 const auto m_i = lumped_mass_matrix->local_element(i);
633 const auto rho_i = density_->local_element(i);
634 const auto e_i = src.local_element(i);
635 dst.local_element(i) = m_i * rho_i * e_i;
636 }
637
638 /* Apply action of diffusion operator \sum_j beta_ij e_j: */
639
640 const auto integrator = [this](const auto &data,
641 auto &dst,
642 const auto &src,
643 const auto range) {
644 dealii::FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(
645 data);
646
647 for (unsigned int cell = range.first; cell < range.second; ++cell) {
648 energy.reinit(cell);
649 energy.read_dof_values(src);
650 apply_local_operator(energy);
651 energy.distribute_local_to_global(dst);
652 }
653 };
654
655 matrix_free_->template cell_loop<vector_type, vector_type>(
656 integrator, dst, src, /* zero destination */ false);
657
658 /* Fix up constrained degrees of freedom: */
659
660 const auto &boundary_map =
661 (level_ == dealii::numbers::invalid_unsigned_int)
662 ? offline_data_->boundary_map()
663 : offline_data_->level_boundary_map()[level_];
664
665 for (auto entry : boundary_map) {
666 const auto i = std::get<0>(entry);
667 if (i >= n_owned)
668 continue;
669
670 const auto id = std::get<4>(entry);
671 if (id == Boundary::dirichlet)
672 dst.local_element(i) = src.local_element(i);
673 }
674 }
675
677 std::shared_ptr<dealii::DiagonalMatrix<vector_type>> &matrix) const
678 {
679 Assert(level_ != dealii::numbers::invalid_unsigned_int,
680 dealii::ExcNotImplemented());
681 matrix = std::make_shared<dealii::DiagonalMatrix<vector_type>>();
682 vector_type &vector = matrix->get_vector();
683 matrix_free_->initialize_dof_vector(vector);
684
685 const vector_type &lumped_mass_matrix =
686 offline_data_->level_lumped_mass_matrix()[level_];
687
688 unsigned int dummy = 0;
689 matrix_free_->template cell_loop<vector_type, unsigned int>(
690 [this](
691 const auto &data, auto &dst, const auto &, const auto range) {
692 dealii::FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(
693 data);
694 dealii::FEEvaluation<dim, order_fe, order_quad, 1, Number> writer(
695 data);
696
697 for (unsigned int cell = range.first; cell < range.second;
698 ++cell) {
699 energy.reinit(cell);
700 writer.reinit(cell);
701 for (unsigned int i = 0; i < energy.dofs_per_cell; ++i) {
702 for (unsigned int j = 0; j < energy.dofs_per_cell; ++j)
703 energy.begin_dof_values()[j] =
704 dealii::VectorizedArray<Number>();
705 energy.begin_dof_values()[i] =
706 dealii::make_vectorized_array<Number>(1.);
707 apply_local_operator(energy);
708 writer.begin_dof_values()[i] = energy.begin_dof_values()[i];
709 }
710 writer.distribute_local_to_global(dst);
711 }
712 },
713 vector,
714 dummy,
715 /* zero destination */ true);
716
717 const unsigned int n_owned =
718 lumped_mass_matrix.get_partitioner()->locally_owned_size();
719
721
723 for (unsigned int i = 0; i < n_owned; ++i) {
724 const auto m_i = lumped_mass_matrix.local_element(i);
725 const auto rho_i = density_->local_element(i);
726 vector.local_element(i) =
727 1. / (m_i * rho_i + vector.local_element(i));
728 }
729
731
732 const auto &boundary_map = offline_data_->level_boundary_map()[level_];
733
734 for (auto entry : boundary_map) {
735 const auto i = std::get<0>(entry);
736 if (i >= n_owned)
737 continue;
738
739 const auto id = std::get<4>(entry);
740 if (id == Boundary::dirichlet)
741 vector.local_element(i) = 1.;
742 }
743 }
744
745 private:
746 const OfflineData<dim, Number2> *offline_data_;
747 const dealii::MatrixFree<dim, Number> *matrix_free_;
748 const dealii::LinearAlgebra::distributed::Vector<Number> *density_;
749 Number factor_;
750 unsigned int level_;
751
752 template <typename Evaluator>
753 void apply_local_operator(Evaluator &energy) const
754 {
755 energy.evaluate(dealii::EvaluationFlags::gradients);
756 for (unsigned int q = 0; q < energy.n_q_points; ++q) {
757 energy.submit_gradient(factor_ * energy.get_gradient(q), q);
758 }
759 energy.integrate(dealii::EvaluationFlags::gradients);
760 }
761 };
762
763
767 template <int dim, typename Number>
768 class MGTransferEnergy : public dealii::MGTransferMatrixFree<dim, Number>
769 {
770 public:
771 void build(const dealii::DoFHandler<dim> &dof_handler,
772 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
773 &matrix_free)
774 {
775 dealii::MGTransferMatrixFree<dim, Number>::build(dof_handler);
776 level_matrix_free_ = &matrix_free;
777 }
778
779 template <typename Number2>
781 const dealii::DoFHandler<dim> &dof_handler,
782 dealii::MGLevelObject<
783 dealii::LinearAlgebra::distributed::Vector<Number>> &dst,
784 const dealii::LinearAlgebra::distributed::Vector<Number2> &src) const
785 {
786 if (dst[dst.min_level()].size() == 0)
787 for (unsigned int l = dst.min_level(); l <= dst.max_level(); ++l)
788 (*level_matrix_free_)[l].initialize_dof_vector(dst[l]);
789 dealii::MGTransferMatrixFree<dim, Number>::copy_to_mg(
790 dof_handler, dst, src);
791 }
792
793 private:
794 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
795 *level_matrix_free_;
796 };
797
798 } // namespace NavierStokes
799} /* namespace ryujin */
800
801#undef locally_owned_size
void vmult(vector_type &dst, const vector_type &src) const
void vmult(block_vector_type &dst, const block_vector_type &src) const
dealii::LinearAlgebra::distributed::Vector< Number > vector_type
void reinit(const vector_type &lumped_mass_matrix, const vector_type &density, const dealii::AffineConstraints< Number > &affine_constraints)
dealii::LinearAlgebra::distributed::BlockVector< Number > block_vector_type
dealii::LinearAlgebra::distributed::Vector< Number > vector_type
void compute_diagonal(std::shared_ptr< dealii::DiagonalMatrix< vector_type > > &matrix) const
void Tvmult(vector_type &dst, const vector_type &src) const
void vmult(vector_type &dst, const vector_type &src) const
Number el(const unsigned int, const unsigned int) const
dealii::types::global_dof_index m() const
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 build(const dealii::DoFHandler< dim > &dof_handler, const dealii::MGLevelObject< dealii::MatrixFree< dim, Number > > &matrix_free)
void copy_to_mg(const dealii::DoFHandler< dim > &dof_handler, dealii::MGLevelObject< dealii::LinearAlgebra::distributed::Vector< Number > > &dst, const dealii::LinearAlgebra::distributed::Vector< Number2 > &src) const
void prolongate(const unsigned int to_level, vector_type &dst, const vector_type &src) const override
void copy_to_mg(const dealii::DoFHandler< dim > &dof_handler, dealii::MGLevelObject< vector_type > &dst, const dealii::LinearAlgebra::distributed::BlockVector< Number2 > &src) const
void restrict_and_add(const unsigned int to_level, vector_type &dst, const vector_type &src) const override
dealii::LinearAlgebra::distributed::Vector< Number > scalar_type
dealii::LinearAlgebra::distributed::BlockVector< Number > vector_type
void interpolate_to_mg(const dealii::DoFHandler< dim > &dof_handler, dealii::MGLevelObject< scalar_type > &dst, const dealii::LinearAlgebra::distributed::Vector< Number2 > &src) const
void build(const dealii::DoFHandler< dim > &dof_handler, const dealii::MGConstrainedDoFs &mg_constrained_dofs, const dealii::MGLevelObject< dealii::MatrixFree< dim, Number > > &matrix_free)
void copy_from_mg(const dealii::DoFHandler< dim > &dof_handler, dealii::LinearAlgebra::distributed::BlockVector< Number2 > &dst, const dealii::MGLevelObject< vector_type > &src) const
void vmult(block_vector_type &dst, const block_vector_type &src) const
dealii::LinearAlgebra::distributed::BlockVector< Number > block_vector_type
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)
dealii::LinearAlgebra::distributed::Vector< Number > vector_type
void Tvmult(block_vector_type &dst, const block_vector_type &src) const
void compute_diagonal(std::shared_ptr< DiagonalMatrix< dim, Number > > &matrix) const
const auto & boundary_map() const
Definition: offline_data.h:196
const auto & level_lumped_mass_matrix() const
Definition: offline_data.h:252
const auto & level_boundary_map() const
Definition: offline_data.h:210
const auto & lumped_mass_matrix() const
Definition: offline_data.h:240
#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
Euler::HyperbolicSystem HyperbolicSystem
Subscriptor EnableObserverPointer