8#include <compile_time_options.h>
15#include "../euler/hyperbolic_system.h"
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>
31 namespace NavierStokes
43 template <
int dim,
typename Number>
50 using vector_type = dealii::LinearAlgebra::distributed::Vector<Number>;
56 dealii::LinearAlgebra::distributed::BlockVector<Number>;
69 const dealii::AffineConstraints<Number> &affine_constraints)
71 diagonal.reinit(density,
true);
73 DEAL_II_OPENMP_SIMD_PRAGMA
74 for (
unsigned int i = 0;
75 i < density.get_partitioner()->locally_owned_size();
77 diagonal.local_element(i) =
79 (density.local_element(i) * lumped_mass_matrix.local_element(i));
86 affine_constraints.set_zero(diagonal);
102 return diagonal_block;
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();
115 dst.local_element(i) =
116 diagonal.local_element(i) * src.local_element(i);
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();
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);
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();
140 dst.block(d).local_element(i) =
141 diagonal_block.block(d).local_element(i) *
142 src.block(d).local_element(i);
158 template <
int dim,
typename Number,
typename Number2>
166 using vector_type = dealii::LinearAlgebra::distributed::Vector<Number>;
168 dealii::LinearAlgebra::distributed::BlockVector<Number>;
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)
180 parabolic_system_ = ¶bolic_system;
181 offline_data_ = &offline_data;
182 matrix_free_ = &matrix_free;
184 theta_x_tau_ = theta_x_tau;
197 using VA = dealii::VectorizedArray<Number>;
198 constexpr auto simd_length = VA::size();
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)
209 Assert(level_ == dealii::numbers::invalid_unsigned_int,
210 dealii::ExcInternalError());
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;
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);
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);
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;
247 const auto integrator = [
this](
const auto &data,
251 dealii::FEEvaluation<dim, order_fe, order_quad, dim, Number> velocity(
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);
262 matrix_free_->template cell_loop<block_vector_type, block_vector_type>(
263 integrator, dst, src,
false);
267 const auto &boundary_map =
268 level_ == dealii::numbers::invalid_unsigned_int
272 for (
auto entry : boundary_map) {
274 const auto i = std::get<0>(entry);
278 const dealii::Tensor<1, dim, Number> normal = std::get<1>(entry);
279 const auto id = std::get<4>(entry);
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);
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;
293 for (
unsigned int d = 0; d < dim; ++d)
294 dst.block(d).local_element(i) = V_i[d];
299 for (
unsigned int d = 0; d < dim; ++d)
300 dst.block(d).local_element(i) = src.block(d).local_element(i);
308 Assert(level_ != dealii::numbers::invalid_unsigned_int,
309 dealii::ExcNotImplemented());
310 matrix = std::make_shared<DiagonalMatrix<dim, Number>>();
313 for (
unsigned int d = 0; d < dim; ++d)
314 matrix_free_->initialize_dof_vector(vector.block(d));
315 vector.collect_sizes();
317 const auto &lumped_mass_matrix =
320 unsigned int dummy = 0;
321 matrix_free_->template cell_loop<block_vector_type, unsigned int>(
323 const auto &data,
auto &dst,
const auto &,
const auto range) {
324 dealii::FEEvaluation<dim, order_fe, order_quad, dim, Number>
326 dealii::FEEvaluation<dim, order_fe, order_quad, dim, Number>
329 for (
unsigned int cell = range.first; cell < range.second;
331 velocity.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];
342 writer.distribute_local_to_global(dst);
349 const unsigned int n_owned =
350 lumped_mass_matrix.get_partitioner()->locally_owned_size();
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));
367 for (
auto entry : boundary_map) {
369 const auto i = std::get<0>(entry);
373 const dealii::Tensor<1, dim, Number> normal = std::get<1>(entry);
374 const auto id = std::get<4>(entry);
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);
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;
387 for (
unsigned int d = 0; d < dim; ++d)
388 vector.block(d).local_element(i) = V_i[d];
393 for (
unsigned int d = 0; d < dim; ++d)
394 vector.block(d).local_element(i) = 1.;
402 const dealii::MatrixFree<dim, Number> *matrix_free_;
407 template <
typename Evaluator>
408 void apply_local_operator(Evaluator &velocity)
const
410 const Number mu = parabolic_system_->mu();
411 const Number lambda = parabolic_system_->lambda();
413 velocity.evaluate(dealii::EvaluationFlags::gradients);
415 for (
unsigned int q = 0; q < velocity.n_q_points; ++q) {
416 if constexpr (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);
424 const auto symmetric_gradient = velocity.get_symmetric_gradient(q);
425 const auto divergence = trace(symmetric_gradient);
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);
434 velocity.integrate(dealii::EvaluationFlags::gradients);
442 template <
int dim,
typename Number>
444 :
public dealii::MGTransferBase<
445 dealii::LinearAlgebra::distributed::BlockVector<Number>>
448 using scalar_type = dealii::LinearAlgebra::distributed::Vector<Number>;
450 dealii::LinearAlgebra::distributed::BlockVector<Number>;
454 void build(
const dealii::DoFHandler<dim> &dof_handler,
455 const dealii::MGConstrainedDoFs &mg_constrained_dofs,
456 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
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();
466 matrix_free[level].initialize_dof_vector(scalar_vector[level]);
473 for (
unsigned int block = 0; block < src.n_blocks(); ++block)
474 transfer_.prolongate(to_level, dst.block(block), src.block(block));
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));
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
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);
498 template <
typename Number2>
501 dealii::MGLevelObject<vector_type> &dst,
502 const dealii::LinearAlgebra::distributed::BlockVector<Number2>
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();
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();
518 dst[level].block(block).copy_locally_owned_data_from(
519 scalar_vector[level]);
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
529 for (
unsigned int block = 0; block < dst.n_blocks(); ++block) {
530 for (
unsigned int level = src.min_level(); level <= src.max_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);
539 dealii::MGTransferMatrixFree<dim, Number> transfer_;
540 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
542 mutable dealii::MGLevelObject<scalar_type> scalar_vector;
549 template <
int dim,
typename Number,
typename Number2>
557 using vector_type = dealii::LinearAlgebra::distributed::Vector<Number>;
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)
568 offline_data_ = &offline_data;
569 matrix_free_ = &matrix_free;
571 factor_ = time_factor;
580 dealii::types::global_dof_index
m()
const
582 return density_->size();
585 Number
el(
const unsigned int,
const unsigned int)
const
587 Assert(
false, dealii::ExcNotImplemented());
595 using VA = dealii::VectorizedArray<Number>;
596 constexpr auto simd_length = VA::size();
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)
607 Assert(level_ == dealii::numbers::invalid_unsigned_int,
608 dealii::ExcInternalError());
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;
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);
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;
640 const auto integrator = [
this](
const auto &data,
644 dealii::FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(
647 for (
unsigned int cell = range.first; cell < range.second; ++cell) {
649 energy.read_dof_values(src);
650 apply_local_operator(energy);
651 energy.distribute_local_to_global(dst);
655 matrix_free_->template cell_loop<vector_type, vector_type>(
656 integrator, dst, src,
false);
660 const auto &boundary_map =
661 (level_ == dealii::numbers::invalid_unsigned_int)
665 for (
auto entry : boundary_map) {
666 const auto i = std::get<0>(entry);
670 const auto id = std::get<4>(entry);
672 dst.local_element(i) = src.local_element(i);
677 std::shared_ptr<dealii::DiagonalMatrix<vector_type>> &matrix)
const
679 Assert(level_ != dealii::numbers::invalid_unsigned_int,
680 dealii::ExcNotImplemented());
681 matrix = std::make_shared<dealii::DiagonalMatrix<vector_type>>();
683 matrix_free_->initialize_dof_vector(vector);
688 unsigned int dummy = 0;
689 matrix_free_->template cell_loop<vector_type, unsigned int>(
691 const auto &data,
auto &dst,
const auto &,
const auto range) {
692 dealii::FEEvaluation<dim, order_fe, order_quad, 1, Number> energy(
694 dealii::FEEvaluation<dim, order_fe, order_quad, 1, Number> writer(
697 for (
unsigned int cell = range.first; cell < range.second;
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];
710 writer.distribute_local_to_global(dst);
717 const unsigned int n_owned =
718 lumped_mass_matrix.get_partitioner()->locally_owned_size();
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));
734 for (
auto entry : boundary_map) {
735 const auto i = std::get<0>(entry);
739 const auto id = std::get<4>(entry);
741 vector.local_element(i) = 1.;
747 const dealii::MatrixFree<dim, Number> *matrix_free_;
748 const dealii::LinearAlgebra::distributed::Vector<Number> *density_;
752 template <
typename Evaluator>
753 void apply_local_operator(Evaluator &energy)
const
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);
759 energy.integrate(dealii::EvaluationFlags::gradients);
767 template <
int dim,
typename Number>
771 void build(
const dealii::DoFHandler<dim> &dof_handler,
772 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
775 dealii::MGTransferMatrixFree<dim, Number>::build(dof_handler);
776 level_matrix_free_ = &matrix_free;
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
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);
794 const dealii::MGLevelObject<dealii::MatrixFree<dim, Number>>
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)
vector_type & get_vector()
dealii::LinearAlgebra::distributed::BlockVector< Number > block_vector_type
block_vector_type & get_block_vector()
static constexpr unsigned int order_quad
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
static constexpr unsigned int order_fe
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
MGTransferVelocity()=default
void vmult(block_vector_type &dst, const block_vector_type &src) const
static constexpr unsigned int order_quad
dealii::LinearAlgebra::distributed::BlockVector< Number > block_vector_type
void initialize(const ParabolicSystem ¶bolic_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
static constexpr unsigned int order_fe
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
const auto & level_lumped_mass_matrix() const
const auto & level_boundary_map() const
const auto & lumped_mass_matrix() const
#define RYUJIN_PARALLEL_REGION_BEGIN
#define RYUJIN_PARALLEL_REGION_END
Euler::HyperbolicSystem HyperbolicSystem
Subscriptor EnableObserverPointer