deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
tensor_product_matrix_creator.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2022 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_tensor_product_matrix_creator_h
16#define dealii_tensor_product_matrix_creator_h
17
18
19#include <deal.II/base/config.h>
20
22
24
25#include <deal.II/fe/fe.h>
26#include <deal.II/fe/fe_dgq.h>
27#include <deal.II/fe/fe_q.h>
28#include <deal.II/fe/fe_tools.h>
31
33#include <deal.II/grid/tria.h>
34
35#include <set>
36
38
39
46{
56
66 template <int dim, typename Number>
67 std::pair<std::array<FullMatrix<Number>, dim>,
68 std::array<FullMatrix<Number>, dim>>
70 const FiniteElement<1> &fe,
71 const Quadrature<1> &quadrature,
72 const ::ndarray<LaplaceBoundaryType, dim, 2> &boundary_ids,
73 const ::ndarray<double, dim, 3> &cell_extent,
74 const unsigned int n_overlap = 1);
75
80 template <int dim, typename Number>
81 std::pair<std::array<FullMatrix<Number>, dim>,
82 std::array<FullMatrix<Number>, dim>>
84 const typename Triangulation<dim>::cell_iterator &cell,
85 const std::set<types::boundary_id> &dirichlet_boundaries,
86 const std::set<types::boundary_id> &neumann_boundaries,
87 const FiniteElement<1> &fe,
88 const Quadrature<1> &quadrature,
89 const ::ndarray<double, dim, 3> &cell_extent,
90 const unsigned int n_overlap = 1);
91
92
93
116 template <typename Number = double>
119 const FiniteElement<1> &fe,
120 const Number &h,
121 const std::pair<bool, bool> include_endpoints = {true, true},
122 std::vector<unsigned int> numbering = std::vector<unsigned int>());
123
124
125
144 template <typename Number = double>
145 FullMatrix<Number>
147 const FiniteElement<1> &fe,
148 const Number &h,
149 const std::pair<bool, bool> include_endpoints = {true, true},
150 std::vector<unsigned int> numbering = std::vector<unsigned int>());
151
152
184 template <typename Number = double>
185 FullMatrix<Number>
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});
191
192
209 template <typename Number = double>
210 FullMatrix<Number>
212 const FiniteElement<1> &fe,
213 const Number h,
214 std::vector<Number> coefficients = std::vector<Number>());
215
216
217
235 template <typename Number = double>
236 FullMatrix<Number>
238 const std::vector<Polynomials::Polynomial<double>>
239 &polynomial_basis_derivative,
240 const unsigned int overlap = 1);
241
242} // namespace TensorProductMatrixCreator
243
244
245
246/*----------------------- Inline functions ----------------------------------*/
247
248
250{
251 namespace internal
252 {
253 template <typename Number>
254 void
255 clear_row_and_column(const unsigned int n_dofs_1D_with_overlap,
256 const unsigned int n,
257 FullMatrix<Number> &matrix)
258 {
259 for (unsigned int i = 0; i < n_dofs_1D_with_overlap; ++i)
260 {
261 matrix[i][n] = 0.0;
262 matrix[n][i] = 0.0;
263 }
264 }
265
266 template <typename Number>
267 std::tuple<FullMatrix<Number>, FullMatrix<Number>, bool>
269 const FiniteElement<1> &fe,
270 const Quadrature<1> &quadrature)
271 {
272 Triangulation<1> tria;
274
275 DoFHandler<1> dof_handler(tria);
276 dof_handler.distribute_dofs(fe);
277
279
280 const unsigned int n_dofs_1D = fe.n_dofs_per_cell();
281
282 FullMatrix<Number> mass_matrix_reference(n_dofs_1D, n_dofs_1D);
283 FullMatrix<Number> derivative_matrix_reference(n_dofs_1D, n_dofs_1D);
284
285 FEValues<1> fe_values(mapping,
286 fe,
287 quadrature,
290
291 fe_values.reinit(tria.begin());
292
293 const auto lexicographic_to_hierarchic_numbering =
296 fe.tensor_degree()));
297
298 for (const unsigned int q_index : fe_values.quadrature_point_indices())
299 for (const unsigned int i : fe_values.dof_indices())
300 for (const unsigned int j : fe_values.dof_indices())
301 {
302 mass_matrix_reference(i, j) +=
303 (fe_values.shape_value(lexicographic_to_hierarchic_numbering[i],
304 q_index) *
305 fe_values.shape_value(lexicographic_to_hierarchic_numbering[j],
306 q_index) *
307 fe_values.JxW(q_index));
308
309 derivative_matrix_reference(i, j) +=
310 (fe_values.shape_grad(lexicographic_to_hierarchic_numbering[i],
311 q_index) *
312 fe_values.shape_grad(lexicographic_to_hierarchic_numbering[j],
313 q_index) *
314 fe_values.JxW(q_index));
315 }
316
317 return std::tuple<FullMatrix<Number>, FullMatrix<Number>, bool>{
318 mass_matrix_reference, derivative_matrix_reference, false};
319 }
320 } // namespace internal
321
322
323
324 template <int dim, typename Number>
325 std::pair<std::array<FullMatrix<Number>, dim>,
326 std::array<FullMatrix<Number>, dim>>
328 const FiniteElement<1> &fe,
329 const Quadrature<1> &quadrature,
330 const ::ndarray<LaplaceBoundaryType, dim, 2> &boundary_ids,
331 const ::ndarray<double, dim, 3> &cell_extent,
332 const unsigned int n_overlap)
333 {
334 // 1) create element mass and siffness matrix (without overlap)
335 const auto create_reference_mass_and_stiffness_matrices =
337 fe, quadrature);
338
339 const auto &M_ref =
340 std::get<0>(create_reference_mass_and_stiffness_matrices);
341 const auto &K_ref =
342 std::get<1>(create_reference_mass_and_stiffness_matrices);
343 const auto &is_dg =
344 std::get<2>(create_reference_mass_and_stiffness_matrices);
345
346 AssertIndexRange(n_overlap, M_ref.n());
347 AssertIndexRange(0, n_overlap);
348 AssertThrow(is_dg == false, ExcNotImplemented());
349
350 // 2) loop over all dimensions and create 1d mass and stiffness
351 // matrices so that boundary conditions and overlap are considered
352
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;
355
356 std::array<FullMatrix<Number>, dim> Ms;
357 std::array<FullMatrix<Number>, dim> Ks;
358
359 for (unsigned int d = 0; d < dim; ++d)
360 {
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);
363
364 // inner cell
365 for (unsigned int i = 0; i < n_dofs_1D; ++i)
366 for (unsigned int j = 0; j < n_dofs_1D; ++j)
367 {
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];
372 }
373
374 // left neighbor or left boundary
375 if (boundary_ids[d][0] == LaplaceBoundaryType::internal_boundary)
376 {
377 // left neighbor
378 Assert(cell_extent[d][0] > 0.0, ExcInternalError());
379
380 for (unsigned int i = 0; i < n_overlap; ++i)
381 for (unsigned int j = 0; j < n_overlap; ++j)
382 {
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];
387 }
388 }
389 else
390 {
391 if (boundary_ids[d][0] == LaplaceBoundaryType::dirichlet)
392 {
393 // left DBC
394 const unsigned i0 = n_overlap - 1;
395 internal::clear_row_and_column(n_dofs_1D_with_overlap,
396 i0,
397 Ms[d]);
398 internal::clear_row_and_column(n_dofs_1D_with_overlap,
399 i0,
400 Ks[d]);
401 }
402 else if (boundary_ids[d][0] == LaplaceBoundaryType::neumann)
403 {
404 // left NBC -> nothing to do
405 }
406 else
407 {
409 }
410 }
411
412 // right neighbor or right boundary
413 if (boundary_ids[d][1] == LaplaceBoundaryType::internal_boundary)
414 {
415 Assert(cell_extent[d][2] > 0.0, ExcInternalError());
416
417 for (unsigned int i = 0; i < n_overlap; ++i)
418 for (unsigned int j = 0; j < n_overlap; ++j)
419 {
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];
424 }
425 }
426 else
427 {
428 if (boundary_ids[d][1] == LaplaceBoundaryType::dirichlet)
429 {
430 // right DBC
431 const unsigned i0 = n_overlap + n_dofs_1D - 2;
432 internal::clear_row_and_column(n_dofs_1D_with_overlap,
433 i0,
434 Ms[d]);
435 internal::clear_row_and_column(n_dofs_1D_with_overlap,
436 i0,
437 Ks[d]);
438 }
439 else if (boundary_ids[d][1] == LaplaceBoundaryType::neumann)
440 {
441 // right NBC -> nothing to do
442 }
443 else
444 {
446 }
447 }
448 }
449
450 return {Ms, Ks};
451 }
452
453
454 template <int dim, typename Number>
455 std::pair<std::array<FullMatrix<Number>, dim>,
456 std::array<FullMatrix<Number>, dim>>
458 const typename Triangulation<dim>::cell_iterator &cell,
459 const std::set<types::boundary_id> &dirichlet_boundaries,
460 const std::set<types::boundary_id> &neumann_boundaries,
461 const FiniteElement<1> &fe,
462 const Quadrature<1> &quadrature,
463 const ::ndarray<double, dim, 3> &cell_extent,
464 const unsigned int n_overlap)
465 {
467
468 for (unsigned int d = 0; d < dim; ++d)
469 {
470 // left neighbor or left boundary
471 if ((cell->at_boundary(2 * d) == false) ||
472 cell->has_periodic_neighbor(2 * d))
473 {
474 // left neighbor
475 Assert(cell_extent[d][0] > 0.0, ExcInternalError());
476
477 boundary_ids[d][0] = LaplaceBoundaryType::internal_boundary;
478 }
479 else
480 {
481 const auto bid = cell->face(2 * d)->boundary_id();
482 if (dirichlet_boundaries.find(bid) !=
483 dirichlet_boundaries.end() /*DBC*/)
484 {
485 // left DBC
486 boundary_ids[d][0] = LaplaceBoundaryType::dirichlet;
487 }
488 else if (neumann_boundaries.find(bid) !=
489 neumann_boundaries.end() /*NBC*/)
490 {
491 // left NBC
492 boundary_ids[d][0] = LaplaceBoundaryType::neumann;
493 }
494 else
495 {
497 }
498 }
499
500 // right neighbor or right boundary
501 if ((cell->at_boundary(2 * d + 1) == false) ||
502 cell->has_periodic_neighbor(2 * d + 1))
503 {
504 Assert(cell_extent[d][2] > 0.0, ExcInternalError());
505
506 boundary_ids[d][1] = LaplaceBoundaryType::internal_boundary;
507 }
508 else
509 {
510 const auto bid = cell->face(2 * d + 1)->boundary_id();
511 if (dirichlet_boundaries.find(bid) !=
512 dirichlet_boundaries.end() /*DBC*/)
513 {
514 // right DBC
515 boundary_ids[d][1] = LaplaceBoundaryType::dirichlet;
516 }
517 else if (neumann_boundaries.find(bid) !=
518 neumann_boundaries.end() /*NBC*/)
519 {
520 // right NBC
521 boundary_ids[d][1] = LaplaceBoundaryType::neumann;
522 }
523 else
524 {
526 }
527 }
528 }
529
531 fe, quadrature, boundary_ids, cell_extent, n_overlap);
532 }
533
534 template <typename Number>
537 const Number &h,
538 const std::pair<bool, bool> include_endpoints,
539 std::vector<unsigned int> numbering)
540 {
541 if (dynamic_cast<const FE_DGQ<1> *>(&fe) == nullptr &&
542 numbering.size() == 0)
543 {
544 Assert(
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."));
549 }
550
551 if (numbering.size() == 0)
552 {
553 numbering.resize(fe.dofs_per_cell);
554 std::iota(numbering.begin(), numbering.end(), 0);
555 }
556
557 const unsigned int degree = fe.degree;
558 const unsigned int n_dofs_per_cell = fe.dofs_per_cell;
559 QGauss<1> quadrature(degree + 1);
560
561 FullMatrix<Number> cell_matrix(n_dofs_per_cell, n_dofs_per_cell);
562 cell_matrix = 0;
563
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;
568
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));
576
577 return cell_matrix;
578 }
579
580 template <typename Number>
583 const Number &h,
584 const std::pair<bool, bool> include_endpoints,
585 std::vector<unsigned int> numbering)
586 {
587 if (dynamic_cast<const FE_DGQ<1> *>(&fe) == nullptr &&
588 numbering.size() == 0)
589 {
590 Assert(
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."));
595 }
596
597 if (numbering.size() == 0)
598 {
599 numbering.resize(fe.dofs_per_cell);
600 std::iota(numbering.begin(), numbering.end(), 0);
601 }
602
603 const unsigned int degree = fe.degree;
604 const unsigned int n_dofs_per_cell = fe.dofs_per_cell;
605 const Number &JxW = h;
606 QGauss<1> quadrature(degree + 1);
607
608 FullMatrix<Number> cell_matrix(n_dofs_per_cell, n_dofs_per_cell);
609 cell_matrix = 0;
610
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;
615
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));
623
624 return cell_matrix;
625 }
626
627 template <typename Number>
630 const unsigned int &n_cells,
631 const unsigned int &overlap,
632 const std::pair<bool, bool> include_endpoints)
633 {
634 const unsigned int n_dofs_per_cell = cell_matrix.n();
635
636 Assert(cell_matrix.m() == n_dofs_per_cell,
638 "The provided cell mass matrix must be a square matrix."));
640 n_cells <= 10,
642 "create_1D_discretization_matrix() returns a full matrix and is not meant to be used with a larger number of cells. "));
643 Assert(n_cells > 0,
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."));
647
648 unsigned int n_total_dofs =
649 n_cells * n_dofs_per_cell - overlap * (n_cells - 1);
650
651 if (!include_endpoints.first)
652 n_total_dofs -= 1;
653 if (!include_endpoints.second)
654 n_total_dofs -= 1;
655
656 FullMatrix<Number> result_matrix(n_total_dofs, n_total_dofs);
657 result_matrix = 0;
658
659 const unsigned int left_shift = include_endpoints.first ? 0 : 1;
660
661 for (unsigned int cell = 0; cell < n_cells; ++cell)
662 {
663 const unsigned int dof_shift = cell * overlap + left_shift;
664
665 const unsigned int start_dof =
666 (cell == 0 && !include_endpoints.first) ? 1 : 0;
667
668 const unsigned int end_dof =
669 (cell == n_cells - 1 && !include_endpoints.second) ?
670 n_dofs_per_cell - 1 :
671 n_dofs_per_cell;
672 for (unsigned int i = start_dof; i < end_dof; ++i)
673 for (unsigned int j = start_dof; j < end_dof; ++j)
674 {
675 result_matrix(i + cell * n_dofs_per_cell - dof_shift,
676 j + cell * n_dofs_per_cell - dof_shift) +=
677 cell_matrix(i, j);
678 }
679 }
680 return result_matrix;
681 }
682
683
684
685 template <typename Number>
688 const Number h,
689 std::vector<Number> coefficients)
690 {
691 Assert(dynamic_cast<const FE_Q<1> *>(&fe) != nullptr, ExcNotImplemented());
692 Assert(h > 0, ExcMessage("Provided element size h is negative"));
693
694 const unsigned int degree = fe.degree;
695 Assert(degree > 0,
696 ExcMessage("Provided element degree has to greater than 0"));
697
698
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"));
703
704 if (coefficients.size() == 0)
705 {
706 coefficients.resize(degree);
707
708 double inverse_factorial_square = 1.;
709 coefficients[0] = 1.;
710 for (unsigned int k = 2; k <= degree; ++k)
711 {
712 inverse_factorial_square /= (k * k);
713 coefficients[k - 1] = inverse_factorial_square;
714 }
715 }
716
717 std::vector<std::vector<Polynomials::Polynomial<double>>> polynomial_basis;
718
719 polynomial_basis.resize(degree + 1);
720
721 auto support_points = fe.get_unit_support_points();
722 std::sort(support_points.begin(),
723 support_points.end(),
724 [](const Point<1> &p, const Point<1> &q) -> bool {
725 return p(0) < q(0);
726 });
727
728 polynomial_basis[0] =
730
731 for (unsigned int k = 1; k < degree + 1; ++k)
732 {
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());
737 }
738
739
740 FullMatrix<Number> penalty_matrix =
741 create_1d_ghost_penalty_matrix(polynomial_basis[1]);
742 penalty_matrix *= coefficients[0];
743
744 for (unsigned int k = 2; k < degree + 1; ++k)
745 {
746 FullMatrix<Number> kth_matrix =
747 create_1d_ghost_penalty_matrix(polynomial_basis[k]);
748 penalty_matrix.add(coefficients[k - 1], kth_matrix);
749 }
750
751 penalty_matrix *= (1 / h);
752 return penalty_matrix;
753 }
754
755
756 template <typename Number>
759 const std::vector<Polynomials::Polynomial<double>>
760 &polynomial_basis_derivative,
761 const unsigned int overlap)
762 {
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;
766
767 FullMatrix<Number> penalty_matrix(n_total_dofs, n_total_dofs);
768
769 std::vector<double> values_left(n_dofs_per_cell);
770 std::vector<double> values_right(n_dofs_per_cell);
771
772 for (unsigned int i = 0; i < n_dofs_per_cell; ++i)
773 {
774 values_left[i] = polynomial_basis_derivative[i].value(0);
775 values_right[i] = polynomial_basis_derivative[i].value(1);
776 }
777
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];
781
782
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];
786
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];
790
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];
794
795 return penalty_matrix;
796 }
797} // namespace TensorProductMatrixCreator
798
799
800
802
803#endif
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)
Definition fe_q.h:554
const unsigned int degree
Definition fe_data.h:452
unsigned int n_dofs_per_cell() const
unsigned int tensor_degree() const
const unsigned int dofs_per_cell
Definition fe_data.h:436
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)
Definition point.h:113
cell_iterator begin(const unsigned int level=0) const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
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
Definition tria.h:1557
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_gradients
Shape function gradients.
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
std::vector< unsigned int > hierarchic_to_lexicographic_numbering(unsigned int degree)
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)
Definition utilities.h:1700
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
Definition ndarray.h:107