44 #if defined DEAL_II_WITH_CGAL || defined DEAL_II_WITH_PARMOONOLITH
46 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
62 &immersed_constraints)
67 ExcMessage(
"This function can only work if dim1<=dim0"));
71 ExcMessage(
"The immersed triangulation can only be a "
72 "parallel::shared::triangulation"));
74 const auto &space_fe = space_dh.
get_fe();
75 const auto &immersed_fe = immersed_dh.
get_fe();
78 const unsigned int n_dofs_per_immersed_cell =
79 immersed_fe.n_dofs_per_cell();
81 const unsigned int n_space_fe_components = space_fe.n_components();
82 const unsigned int n_immersed_fe_components = immersed_fe.n_components();
85 n_dofs_per_immersed_cell);
87 std::vector<types::global_dof_index> local_space_dof_indices(
88 n_dofs_per_space_cell);
89 std::vector<types::global_dof_index> local_immersed_dof_indices(
90 n_dofs_per_immersed_cell);
96 (immersed_comps.
size() == 0 ?
103 std::vector<unsigned int> space_gtl(n_space_fe_components,
105 std::vector<unsigned int> immersed_gtl(n_immersed_fe_components,
107 for (
unsigned int i = 0, j = 0; i < n_space_fe_components; ++i)
114 for (
unsigned int i = 0, j = 0; i < n_immersed_fe_components; ++i)
117 immersed_gtl[i] = j++;
122 Table<2, bool> dof_mask(n_dofs_per_space_cell, n_dofs_per_immersed_cell);
123 dof_mask.fill(
false);
125 for (
unsigned int i = 0; i < n_dofs_per_space_cell; ++i)
127 const auto comp_i = space_fe.system_to_component_index(i).first;
130 for (
unsigned int j = 0; j < n_dofs_per_immersed_cell; ++j)
133 immersed_fe.system_to_component_index(j).first;
134 if (immersed_gtl[comp_j] == space_gtl[comp_i])
136 dof_mask(i, j) =
true;
148 for (
const auto &infos : cells_and_quads)
150 const auto &[first_cell, second_cell, quad_formula] = infos;
152 local_cell_matrix =
typename Matrix::value_type();
154 const unsigned int n_quad_pts = quad_formula.size();
155 const auto &real_qpts = quad_formula.get_points();
156 std::vector<Point<dim0>> ref_pts_space(n_quad_pts);
157 std::vector<Point<dim1>> ref_pts_immersed(n_quad_pts);
165 const auto &JxW = quad_formula.get_weights();
166 for (
unsigned int q = 0; q < n_quad_pts; ++q)
168 for (
unsigned int i = 0; i < n_dofs_per_space_cell; ++i)
170 const unsigned int comp_i =
174 for (
unsigned int j = 0; j < n_dofs_per_immersed_cell;
177 const unsigned int comp_j =
181 if (space_gtl[comp_i] == immersed_gtl[comp_j])
183 local_cell_matrix(i, j) +=
184 space_fe.shape_value(i, ref_pts_space[q]) *
185 immersed_fe.shape_value(j,
186 ref_pts_immersed[q]) *
194 *first_cell, &space_dh);
196 *second_cell, &immersed_dh);
197 space_cell_dh->get_dof_indices(local_space_dof_indices);
198 immersed_cell_dh->get_dof_indices(local_immersed_dof_indices);
234 local_space_dof_indices,
235 immersed_constraints,
236 local_immersed_dof_indices,
249 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
268 "This function needs CGAL or PARMOONOLITH to be installed, "
269 "but cmake could not either."));
283 std::tuple<
typename dealii::Triangulation<1, 1>::cell_iterator,
284 typename dealii::Triangulation<1, 1>::cell_iterator,
285 dealii::Quadrature<1>>> &,
286 dealii::SparseMatrix<double> &,
299 std::tuple<
typename dealii::Triangulation<2, 2>::cell_iterator,
300 typename dealii::Triangulation<1, 2>::cell_iterator,
301 dealii::Quadrature<2>>> &,
302 dealii::SparseMatrix<double> &,
315 std::tuple<
typename dealii::Triangulation<2, 2>::cell_iterator,
316 typename dealii::Triangulation<2, 2>::cell_iterator,
317 dealii::Quadrature<2>>> &,
318 dealii::SparseMatrix<double> &,
333 std::tuple<
typename dealii::Triangulation<3, 3>::cell_iterator,
334 typename dealii::Triangulation<1, 3>::cell_iterator,
335 dealii::Quadrature<3>>> &,
336 dealii::SparseMatrix<double> &,
349 std::tuple<
typename dealii::Triangulation<3, 3>::cell_iterator,
350 typename dealii::Triangulation<2, 3>::cell_iterator,
351 dealii::Quadrature<3>>> &,
352 dealii::SparseMatrix<double> &,
380 dealii::TrilinosWrappers::SparseMatrix>(
381 dealii::DoFHandler<2, 2>
const &,
382 dealii::DoFHandler<1, 2>
const &,
384 std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
385 dealii::Triangulation<1, 2>::cell_iterator,
386 dealii::Quadrature<2>>,
387 std::allocator<std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
388 dealii::Triangulation<1, 2>::cell_iterator,
389 dealii::Quadrature<2>>>>
const &,
390 dealii::TrilinosWrappers::SparseMatrix &,
391 dealii::AffineConstraints<
392 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &,
393 dealii::ComponentMask
const &,
394 dealii::ComponentMask
const &,
395 dealii::Mapping<2, 2>
const &,
396 dealii::Mapping<1, 2>
const &,
397 dealii::AffineConstraints<
398 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &);
405 dealii::TrilinosWrappers::SparseMatrix>(
406 dealii::DoFHandler<2, 2>
const &,
407 dealii::DoFHandler<2, 2>
const &,
409 std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
410 dealii::Triangulation<2, 2>::cell_iterator,
411 dealii::Quadrature<2>>,
412 std::allocator<std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
413 dealii::Triangulation<2, 2>::cell_iterator,
414 dealii::Quadrature<2>>>>
const &,
415 dealii::TrilinosWrappers::SparseMatrix &,
416 dealii::AffineConstraints<
417 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &,
418 dealii::ComponentMask
const &,
419 dealii::ComponentMask
const &,
420 dealii::Mapping<2, 2>
const &,
421 dealii::Mapping<2, 2>
const &,
422 dealii::AffineConstraints<
423 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &);
430 dealii::TrilinosWrappers::SparseMatrix>(
431 dealii::DoFHandler<3, 3>
const &,
432 dealii::DoFHandler<2, 3>
const &,
434 std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
435 dealii::Triangulation<2, 3>::cell_iterator,
436 dealii::Quadrature<3>>,
437 std::allocator<std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
438 dealii::Triangulation<2, 3>::cell_iterator,
439 dealii::Quadrature<3>>>>
const &,
440 dealii::TrilinosWrappers::SparseMatrix &,
441 dealii::AffineConstraints<
442 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &,
443 dealii::ComponentMask
const &,
444 dealii::ComponentMask
const &,
445 dealii::Mapping<3, 3>
const &,
446 dealii::Mapping<2, 3>
const &,
447 dealii::AffineConstraints<
448 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &);
455 dealii::TrilinosWrappers::SparseMatrix>(
456 dealii::DoFHandler<3, 3>
const &,
457 dealii::DoFHandler<3, 3>
const &,
459 std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
460 dealii::Triangulation<3, 3>::cell_iterator,
461 dealii::Quadrature<3>>,
462 std::allocator<std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
463 dealii::Triangulation<3, 3>::cell_iterator,
464 dealii::Quadrature<3>>>>
const &,
465 dealii::TrilinosWrappers::SparseMatrix &,
466 dealii::AffineConstraints<
467 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &,
468 dealii::ComponentMask
const &,
469 dealii::ComponentMask
const &,
470 dealii::Mapping<3, 3>
const &,
471 dealii::Mapping<3, 3>
const &,
472 dealii::AffineConstraints<
473 dealii::TrilinosWrappers::SparseMatrix::value_type>
const &);
480 dealii::PETScWrappers::MPI::SparseMatrix>(
481 dealii::DoFHandler<2, 2>
const &,
482 dealii::DoFHandler<1, 2>
const &,
484 std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
485 dealii::Triangulation<1, 2>::cell_iterator,
486 dealii::Quadrature<2>>,
487 std::allocator<std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
488 dealii::Triangulation<1, 2>::cell_iterator,
489 dealii::Quadrature<2>>>>
const &,
490 dealii::PETScWrappers::MPI::SparseMatrix &,
491 dealii::AffineConstraints<
492 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &,
493 dealii::ComponentMask
const &,
494 dealii::ComponentMask
const &,
495 dealii::Mapping<2, 2>
const &,
496 dealii::Mapping<1, 2>
const &,
497 dealii::AffineConstraints<
498 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &);
505 dealii::PETScWrappers::MPI::SparseMatrix>(
506 dealii::DoFHandler<2, 2>
const &,
507 dealii::DoFHandler<2, 2>
const &,
509 std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
510 dealii::Triangulation<2, 2>::cell_iterator,
511 dealii::Quadrature<2>>,
512 std::allocator<std::tuple<dealii::Triangulation<2, 2>::cell_iterator,
513 dealii::Triangulation<2, 2>::cell_iterator,
514 dealii::Quadrature<2>>>>
const &,
515 dealii::PETScWrappers::MPI::SparseMatrix &,
516 dealii::AffineConstraints<
517 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &,
518 dealii::ComponentMask
const &,
519 dealii::ComponentMask
const &,
520 dealii::Mapping<2, 2>
const &,
521 dealii::Mapping<2, 2>
const &,
522 dealii::AffineConstraints<
523 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &);
530 dealii::PETScWrappers::MPI::SparseMatrix>(
531 dealii::DoFHandler<3, 3>
const &,
532 dealii::DoFHandler<2, 3>
const &,
534 std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
535 dealii::Triangulation<2, 3>::cell_iterator,
536 dealii::Quadrature<3>>,
537 std::allocator<std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
538 dealii::Triangulation<2, 3>::cell_iterator,
539 dealii::Quadrature<3>>>>
const &,
540 dealii::PETScWrappers::MPI::SparseMatrix &,
541 dealii::AffineConstraints<
542 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &,
543 dealii::ComponentMask
const &,
544 dealii::ComponentMask
const &,
545 dealii::Mapping<3, 3>
const &,
546 dealii::Mapping<2, 3>
const &,
547 dealii::AffineConstraints<
548 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &);
555 dealii::PETScWrappers::MPI::SparseMatrix>(
556 dealii::DoFHandler<3, 3>
const &,
557 dealii::DoFHandler<3, 3>
const &,
559 std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
560 dealii::Triangulation<3, 3>::cell_iterator,
561 dealii::Quadrature<3>>,
562 std::allocator<std::tuple<dealii::Triangulation<3, 3>::cell_iterator,
563 dealii::Triangulation<3, 3>::cell_iterator,
564 dealii::Quadrature<3>>>>
const &,
565 dealii::PETScWrappers::MPI::SparseMatrix &,
566 dealii::AffineConstraints<
567 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &,
568 dealii::ComponentMask
const &,
569 dealii::ComponentMask
const &,
570 dealii::Mapping<3, 3>
const &,
571 dealii::Mapping<3, 3>
const &,
572 dealii::AffineConstraints<
573 dealii::PETScWrappers::MPI::SparseMatrix::value_type>
const &);
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
unsigned int size() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
types::global_dof_index n_dofs() const
unsigned int n_dofs_per_cell() const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
virtual void transform_points_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const ArrayView< const Point< spacedim > > &real_points, const ArrayView< Point< dim > > &unit_points) const
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
void assemble_coupling_mass_matrix_with_exact_intersections(const DoFHandler< dim0, spacedim > &, const DoFHandler< dim1, spacedim > &, const std::vector< std::tuple< typename Triangulation< dim0, spacedim >::cell_iterator, typename Triangulation< dim1, spacedim >::cell_iterator, Quadrature< spacedim >>> &, Matrix &matrix, const AffineConstraints< typename Matrix::value_type > &, const ComponentMask &, const ComponentMask &, const Mapping< dim0, spacedim > &, const Mapping< dim1, spacedim > &, const AffineConstraints< typename Matrix::value_type > &)
Create the coupling mass matrix for non-matching, overlapping grids in an "exact" way,...
static const unsigned int invalid_unsigned_int