15#ifndef dealii_tensor_product_matrix_creator_h
16#define dealii_tensor_product_matrix_creator_h
66 template <
int dim,
typename Number>
67 std::pair<std::array<FullMatrix<Number>, dim>,
68 std::array<FullMatrix<Number>, dim>>
72 const ::ndarray<LaplaceBoundaryType, dim, 2> &boundary_ids,
73 const ::ndarray<double, dim, 3> &cell_extent,
74 const unsigned int n_overlap = 1);
80 template <
int dim,
typename Number>
81 std::pair<std::array<FullMatrix<Number>, dim>,
82 std::array<FullMatrix<Number>, dim>>
85 const std::set<types::boundary_id> &dirichlet_boundaries,
86 const std::set<types::boundary_id> &neumann_boundaries,
89 const ::ndarray<double, dim, 3> &cell_extent,
90 const unsigned int n_overlap = 1);
116 template <
typename Number =
double>
121 const std::pair<bool, bool> include_endpoints = {
true,
true},
122 std::vector<unsigned int> numbering = std::vector<unsigned int>());
144 template <
typename Number =
double>
147 const FiniteElement<1> &fe,
149 const std::pair<bool, bool> include_endpoints = {
true,
true},
150 std::vector<unsigned int> numbering = std::vector<unsigned int>());
184 template <
typename Number =
double>
187 FullMatrix<Number> &cell_matrix,
188 const unsigned int &n_cells,
189 const unsigned int &overlap,
190 const std::pair<bool, bool> include_endpoints = {
true,
true});
209 template <
typename Number =
double>
212 const FiniteElement<1> &fe,
214 std::vector<Number> coefficients = std::vector<Number>());
235 template <
typename Number =
double>
238 const std::vector<Polynomials::Polynomial<double>>
239 &polynomial_basis_derivative,
240 const unsigned int overlap = 1);
253 template <
typename Number>
256 const unsigned int n,
259 for (
unsigned int i = 0; i < n_dofs_1D_with_overlap; ++i)
266 template <
typename Number>
293 const auto lexicographic_to_hierarchic_numbering =
299 for (
const unsigned int i : fe_values.
dof_indices())
300 for (
const unsigned int j : fe_values.
dof_indices())
302 mass_matrix_reference(i, j) +=
303 (fe_values.
shape_value(lexicographic_to_hierarchic_numbering[i],
305 fe_values.
shape_value(lexicographic_to_hierarchic_numbering[j],
307 fe_values.
JxW(q_index));
309 derivative_matrix_reference(i, j) +=
310 (fe_values.
shape_grad(lexicographic_to_hierarchic_numbering[i],
312 fe_values.
shape_grad(lexicographic_to_hierarchic_numbering[j],
314 fe_values.
JxW(q_index));
318 mass_matrix_reference, derivative_matrix_reference,
false};
324 template <
int dim,
typename Number>
325 std::pair<std::array<FullMatrix<Number>, dim>,
326 std::array<FullMatrix<Number>, dim>>
330 const ::ndarray<LaplaceBoundaryType, dim, 2> &boundary_ids,
331 const ::ndarray<double, dim, 3> &cell_extent,
332 const unsigned int n_overlap)
335 const auto create_reference_mass_and_stiffness_matrices =
340 std::get<0>(create_reference_mass_and_stiffness_matrices);
342 std::get<1>(create_reference_mass_and_stiffness_matrices);
344 std::get<2>(create_reference_mass_and_stiffness_matrices);
353 const unsigned int n_dofs_1D = M_ref.n();
354 const unsigned int n_dofs_1D_with_overlap = M_ref.n() - 2 + 2 * n_overlap;
356 std::array<FullMatrix<Number>, dim> Ms;
357 std::array<FullMatrix<Number>, dim> Ks;
359 for (
unsigned int d = 0; d < dim; ++d)
361 Ms[d].reinit(n_dofs_1D_with_overlap, n_dofs_1D_with_overlap);
362 Ks[d].reinit(n_dofs_1D_with_overlap, n_dofs_1D_with_overlap);
365 for (
unsigned int i = 0; i < n_dofs_1D; ++i)
366 for (
unsigned int j = 0; j < n_dofs_1D; ++j)
368 const unsigned int i0 = i + n_overlap - 1;
369 const unsigned int j0 = j + n_overlap - 1;
370 Ms[d][i0][j0] = M_ref[i][j] * cell_extent[d][1];
371 Ks[d][i0][j0] = K_ref[i][j] / cell_extent[d][1];
380 for (
unsigned int i = 0; i < n_overlap; ++i)
381 for (
unsigned int j = 0; j < n_overlap; ++j)
383 const unsigned int i0 = n_dofs_1D - n_overlap + i;
384 const unsigned int j0 = n_dofs_1D - n_overlap + j;
385 Ms[d][i][j] += M_ref[i0][j0] * cell_extent[d][0];
386 Ks[d][i][j] += K_ref[i0][j0] / cell_extent[d][0];
394 const unsigned i0 = n_overlap - 1;
417 for (
unsigned int i = 0; i < n_overlap; ++i)
418 for (
unsigned int j = 0; j < n_overlap; ++j)
420 const unsigned int i0 = n_overlap + n_dofs_1D + i - 2;
421 const unsigned int j0 = n_overlap + n_dofs_1D + j - 2;
422 Ms[d][i0][j0] += M_ref[i][j] * cell_extent[d][2];
423 Ks[d][i0][j0] += K_ref[i][j] / cell_extent[d][2];
431 const unsigned i0 = n_overlap + n_dofs_1D - 2;
454 template <
int dim,
typename Number>
455 std::pair<std::array<FullMatrix<Number>, dim>,
456 std::array<FullMatrix<Number>, dim>>
459 const std::set<types::boundary_id> &dirichlet_boundaries,
460 const std::set<types::boundary_id> &neumann_boundaries,
463 const ::ndarray<double, dim, 3> &cell_extent,
464 const unsigned int n_overlap)
468 for (
unsigned int d = 0; d < dim; ++d)
481 const auto bid = cell->
face(2 * d)->boundary_id();
482 if (dirichlet_boundaries.find(bid) !=
483 dirichlet_boundaries.end() )
488 else if (neumann_boundaries.find(bid) !=
489 neumann_boundaries.end() )
510 const auto bid = cell->
face(2 * d + 1)->boundary_id();
511 if (dirichlet_boundaries.find(bid) !=
512 dirichlet_boundaries.end() )
517 else if (neumann_boundaries.find(bid) !=
518 neumann_boundaries.end() )
531 fe, quadrature, boundary_ids, cell_extent, n_overlap);
534 template <
typename Number>
538 const std::pair<bool, bool> include_endpoints,
539 std::vector<unsigned int> numbering)
541 if (
dynamic_cast<const FE_DGQ<1> *
>(&fe) ==
nullptr &&
542 numbering.size() == 0)
545 include_endpoints.first ==
true && include_endpoints.second ==
true,
547 "You tried to generate a 1D mass matrix with excluding boundary "
548 "dofs for a non-DGQ element without providing a numbering."));
551 if (numbering.size() == 0)
554 std::iota(numbering.begin(), numbering.end(), 0);
557 const unsigned int degree = fe.
degree;
564 unsigned int start_dof = include_endpoints.first ? 0 : 1;
565 unsigned int end_dof =
566 include_endpoints.second ? n_dofs_per_cell : n_dofs_per_cell - 1;
567 const unsigned int shift = include_endpoints.first ? 0 : 1;
569 for (
unsigned int i = start_dof; i < end_dof; ++i)
570 for (
unsigned int j = start_dof; j < end_dof; ++j)
571 for (
unsigned int q = 0; q < quadrature.size(); ++q)
572 cell_matrix(i - shift, j - shift) +=
573 (fe.
shape_value(numbering[i], quadrature.point(q)) *
574 fe.
shape_value(numbering[j], quadrature.point(q))) *
575 (h * quadrature.weight(q));
580 template <
typename Number>
584 const std::pair<bool, bool> include_endpoints,
585 std::vector<unsigned int> numbering)
587 if (
dynamic_cast<const FE_DGQ<1> *
>(&fe) ==
nullptr &&
588 numbering.size() == 0)
591 include_endpoints.first ==
true && include_endpoints.second ==
true,
593 "You tried to generate a 1D derivative matrix with excluding boundary "
594 "dofs for a non-DGQ element without providing a numbering."));
597 if (numbering.size() == 0)
600 std::iota(numbering.begin(), numbering.end(), 0);
603 const unsigned int degree = fe.
degree;
605 const Number &JxW = h;
611 unsigned int start_dof = include_endpoints.first ? 0 : 1;
612 unsigned int end_dof =
613 include_endpoints.second ? n_dofs_per_cell : n_dofs_per_cell - 1;
614 const unsigned int shift = include_endpoints.first ? 0 : 1;
616 for (
unsigned int i = start_dof; i < end_dof; ++i)
617 for (
unsigned int j = start_dof; j < end_dof; ++j)
618 for (
unsigned int q = 0; q < quadrature.size(); ++q)
619 cell_matrix(i - shift, j - shift) +=
620 (fe.
shape_grad(numbering[i], quadrature.point(q)) / h *
621 fe.
shape_grad(numbering[j], quadrature.point(q))) /
622 h * (h * quadrature.weight(q));
627 template <
typename Number>
630 const unsigned int &n_cells,
631 const unsigned int &overlap,
632 const std::pair<bool, bool> include_endpoints)
634 const unsigned int n_dofs_per_cell = cell_matrix.n();
636 Assert(cell_matrix.m() == n_dofs_per_cell,
638 "The provided cell mass matrix must be a square matrix."));
642 "create_1D_discretization_matrix() returns a full matrix and is not meant to be used with a larger number of cells. "));
644 ExcMessage(
"You are trying to get a mass matrix of zero cells."));
645 Assert(overlap < n_dofs_per_cell,
646 ExcMessage(
"The overlap must be smaller than the number of dofs."));
648 unsigned int n_total_dofs =
649 n_cells * n_dofs_per_cell - overlap * (n_cells - 1);
651 if (!include_endpoints.first)
653 if (!include_endpoints.second)
659 const unsigned int left_shift = include_endpoints.first ? 0 : 1;
661 for (
unsigned int cell = 0; cell < n_cells; ++cell)
663 const unsigned int dof_shift = cell * overlap + left_shift;
665 const unsigned int start_dof =
666 (cell == 0 && !include_endpoints.first) ? 1 : 0;
668 const unsigned int end_dof =
669 (cell == n_cells - 1 && !include_endpoints.second) ?
670 n_dofs_per_cell - 1 :
672 for (
unsigned int i = start_dof; i < end_dof; ++i)
673 for (
unsigned int j = start_dof; j < end_dof; ++j)
675 result_matrix(i + cell * n_dofs_per_cell - dof_shift,
676 j + cell * n_dofs_per_cell - dof_shift) +=
680 return result_matrix;
685 template <
typename Number>
689 std::vector<Number> coefficients)
694 const unsigned int degree = fe.
degree;
696 ExcMessage(
"Provided element degree has to greater than 0"));
699 Assert(coefficients.size() == 0 || coefficients.size() == degree,
701 "Provided coefficients vector has to be empty or the same size "
702 "as the number of dofs"));
704 if (coefficients.size() == 0)
706 coefficients.resize(degree);
708 double inverse_factorial_square = 1.;
709 coefficients[0] = 1.;
710 for (
unsigned int k = 2; k <= degree; ++k)
712 inverse_factorial_square /= (k * k);
713 coefficients[k - 1] = inverse_factorial_square;
717 std::vector<std::vector<Polynomials::Polynomial<double>>> polynomial_basis;
719 polynomial_basis.resize(degree + 1);
722 std::sort(support_points.begin(),
723 support_points.end(),
728 polynomial_basis[0] =
731 for (
unsigned int k = 1; k < degree + 1; ++k)
733 polynomial_basis[k].reserve(degree + 1);
734 for (
unsigned int i = 0; i < degree + 1; ++i)
735 polynomial_basis[k].push_back(
736 polynomial_basis[k - 1][i].derivative());
742 penalty_matrix *= coefficients[0];
744 for (
unsigned int k = 2; k < degree + 1; ++k)
748 penalty_matrix.
add(coefficients[k - 1], kth_matrix);
751 penalty_matrix *= (1 / h);
752 return penalty_matrix;
756 template <
typename Number>
760 &polynomial_basis_derivative,
761 const unsigned int overlap)
763 const unsigned int n_dofs_per_cell = polynomial_basis_derivative.size();
764 const unsigned int n_total_dofs = 2 * n_dofs_per_cell - overlap;
765 const unsigned int shift = n_dofs_per_cell - overlap;
769 std::vector<double> values_left(n_dofs_per_cell);
770 std::vector<double> values_right(n_dofs_per_cell);
772 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
774 values_left[i] = polynomial_basis_derivative[i].value(0);
775 values_right[i] = polynomial_basis_derivative[i].value(1);
778 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
779 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j)
780 penalty_matrix(i, j) += values_right[i] * values_right[j];
783 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
784 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j)
785 penalty_matrix(i + shift, j) -= values_left[i] * values_right[j];
787 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
788 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j)
789 penalty_matrix(i, j + shift) -= values_right[i] * values_left[j];
791 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
792 for (
unsigned int j = 0; j < n_dofs_per_cell; ++j)
793 penalty_matrix(i + shift, j + shift) += values_left[i] * values_left[j];
795 return penalty_matrix;
bool at_boundary(const unsigned int i) const
bool has_periodic_neighbor(const unsigned int i) const
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face(const unsigned int i) const
void distribute_dofs(const FiniteElement< dim, spacedim > &fe)
std_cxx20::ranges::iota_view< unsigned int, unsigned int > dof_indices() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > quadrature_point_indices() const
const Tensor< 1, spacedim > & shape_grad(const unsigned int i, const unsigned int q_point) const
double JxW(const unsigned int q_point) const
const double & shape_value(const unsigned int i, const unsigned int q_point) const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access > > &cell)
const unsigned int degree
unsigned int n_dofs_per_cell() const
unsigned int tensor_degree() const
const unsigned int dofs_per_cell
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
const std::vector< Point< dim > > & get_unit_support_points() const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
void add(const number a, const FullMatrix< number2 > &A)
cell_iterator begin(const unsigned int level=0) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
std::vector< Polynomial< double > > generate_complete_Lagrange_basis(const std::vector< Point< 1 > > &points)
void clear_row_and_column(const unsigned int n_dofs_1D_with_overlap, const unsigned int n, FullMatrix< Number > &matrix)
std::tuple< FullMatrix< Number >, FullMatrix< Number >, bool > create_reference_mass_and_stiffness_matrices(const FiniteElement< 1 > &fe, const Quadrature< 1 > &quadrature)
std::pair< std::array< FullMatrix< Number >, dim >, std::array< FullMatrix< Number >, dim > > create_laplace_tensor_product_matrix(const FiniteElement< 1 > &fe, const Quadrature< 1 > &quadrature, const ::ndarray< LaplaceBoundaryType, dim, 2 > &boundary_ids, const ::ndarray< double, dim, 3 > &cell_extent, const unsigned int n_overlap=1)
FullMatrix< Number > create_1d_cell_mass_matrix(const FiniteElement< 1 > &fe, const Number &h, const std::pair< bool, bool > include_endpoints={true, true}, std::vector< unsigned int > numbering=std::vector< unsigned int >())
FullMatrix< Number > create_1D_discretization_matrix(FullMatrix< Number > &cell_matrix, const unsigned int &n_cells, const unsigned int &overlap, const std::pair< bool, bool > include_endpoints={true, true})
FullMatrix< Number > create_1d_cell_laplace_matrix(const FiniteElement< 1 > &fe, const Number &h, const std::pair< bool, bool > include_endpoints={true, true}, std::vector< unsigned int > numbering=std::vector< unsigned int >())
FullMatrix< Number > create_1d_ghost_penalty_matrix(const FiniteElement< 1 > &fe, const Number h, std::vector< Number > coefficients=std::vector< Number >())
Create a 1D ghost penalty matrix for a given finite element. Ghost penalty is used for stabilization ...
std::vector< Integer > invert_permutation(const std::vector< Integer > &permutation)
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray