63 template <
int dim0,
int dim1,
int spacedim>
64 std::pair<std::vector<Point<spacedim>>, std::vector<unsigned int>>
70 const bool tria_is_parallel)
72 const auto & immersed_fe = immersed_dh.
get_fe();
73 std::vector<Point<spacedim>> points_over_local_cells;
75 std::vector<unsigned int> used_cells_ids;
81 unsigned int cell_id = 0;
84 bool use_cell =
false;
87 const auto bbox = cell->bounding_box();
88 std::vector<std::pair<
93 boost::geometry::index::intersects(bbox),
94 std::back_inserter(out_vals));
99 for (
const auto &bbox_it : out_vals)
100 if (bbox_it.second->is_locally_owned())
103 used_cells_ids.emplace_back(cell_id);
115 const std::vector<Point<spacedim>> &x_points =
116 fe_v.get_quadrature_points();
119 points_over_local_cells.insert(points_over_local_cells.end(),
126 return {std::move(points_over_local_cells), std::move(used_cells_ids)};
136 template <
int dim0,
int dim1,
int spacedim>
137 std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
159 for (
unsigned int i = 0, j = 0; i < gtl0.size(); ++i)
163 for (
unsigned int i = 0, j = 0; i < gtl1.size(); ++i)
199 immersed_constraints);
225 ExcMessage(
"This function can only work if dim1 <= dim0"));
231 const bool tria_is_parallel =
234 const auto &space_fe = space_dh.
get_fe();
235 const auto &immersed_fe = immersed_dh.
get_fe();
238 std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
239 std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
247 (immersed_comps.
size() == 0 ?
255 std::vector<unsigned int> space_gtl(space_fe.n_components(),
257 std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
260 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
264 for (
unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
266 immersed_gtl[i] = j++;
268 const unsigned int n_q_points = quad.
size();
269 const unsigned int n_active_c =
273 cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
275 const auto &points_over_local_cells = std::get<0>(qpoints_cells_data);
276 const auto &used_cells_ids = std::get<1>(qpoints_cells_data);
302 const auto &all_cells = std::get<0>(cpm);
303 const auto &maps = std::get<2>(cpm);
306 std::set<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
307 cell_sets(n_active_c);
309 for (
unsigned int i = 0; i < maps.size(); ++i)
314 unsigned int last_id = std::numeric_limits<unsigned int>::max();
315 unsigned int cell_id;
316 for (
const unsigned int idx : maps[i])
319 if (tria_is_parallel)
320 cell_id = used_cells_ids[idx / n_q_points];
322 cell_id = idx / n_q_points;
324 if (last_id != cell_id)
326 cell_sets[cell_id].insert(all_cells[i]);
338 cell->get_dof_indices(dofs);
341 const auto &cells = cell_sets[i];
343 for (
const auto &cell_c : cells)
349 if (ocell->is_locally_owned())
351 ocell->get_dof_indices(odofs);
357 immersed_constraints,
368 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
393 immersed_constraints);
398 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
415 ExcMessage(
"This function can only work if dim1 <= dim0"));
421 const bool tria_is_parallel =
425 const auto &space_fe = space_dh.
get_fe();
426 const auto &immersed_fe = immersed_dh.
get_fe();
429 std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
430 std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
438 (immersed_comps.
size() == 0 ?
445 std::vector<unsigned int> space_gtl(space_fe.n_components(),
447 std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
450 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
454 for (
unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
456 immersed_gtl[i] = j++;
468 const unsigned int n_q_points = quad.
size();
469 const unsigned int n_active_c =
473 cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
475 const auto &points_over_local_cells = std::get<0>(used_cells_data);
476 const auto &used_cells_ids = std::get<1>(used_cells_data);
481 const auto &all_cells = std::get<0>(cpm);
482 const auto &all_qpoints = std::get<1>(cpm);
483 const auto &all_maps = std::get<2>(cpm);
486 std::vector<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
487 cell_container(n_active_c);
488 std::vector<std::vector<std::vector<Point<dim0>>>> qpoints_container(
490 std::vector<std::vector<std::vector<unsigned int>>> maps_container(
495 for (
unsigned int o = 0; o < all_cells.size(); ++o)
497 for (
unsigned int j = 0; j < all_maps[o].size(); ++j)
502 unsigned int cell_id;
503 if (tria_is_parallel)
504 cell_id = used_cells_ids[all_maps[o][j] / n_q_points];
506 cell_id = all_maps[o][j] / n_q_points;
508 const unsigned int n_pt = all_maps[o][j] % n_q_points;
511 if (cell_container[cell_id].empty())
513 cell_container[cell_id].emplace_back(all_cells[o]);
514 qpoints_container[cell_id].emplace_back(
516 maps_container[cell_id].emplace_back(
517 std::vector<unsigned int>{n_pt});
521 else if (cell_container[cell_id].back() == all_cells[o])
523 qpoints_container[cell_id].back().emplace_back(
525 maps_container[cell_id].back().emplace_back(n_pt);
530 const auto cell_p = std::find(cell_container[cell_id].begin(),
531 cell_container[cell_id].end() - 1,
534 if (cell_p == cell_container[cell_id].end() - 1)
536 cell_container[cell_id].emplace_back(all_cells[o]);
537 qpoints_container[cell_id].emplace_back(
539 maps_container[cell_id].emplace_back(
540 std::vector<unsigned int>{n_pt});
544 const unsigned int pos =
545 cell_p - cell_container[cell_id].begin();
546 qpoints_container[cell_id][pos].emplace_back(
548 maps_container[cell_id][pos].emplace_back(n_pt);
556 endc = immersed_dh.
end();
558 for (
unsigned int j = 0; cell != endc; ++cell, ++j)
562 cell->get_dof_indices(dofs);
565 const auto &cells = cell_container[j];
566 const auto &qpoints = qpoints_container[j];
567 const auto &maps = maps_container[j];
569 for (
unsigned int c = 0; c < cells.size(); ++c)
573 *cells[c], &space_dh);
575 if (ocell->is_locally_owned())
577 const std::vector<Point<dim0>> & qps = qpoints[c];
578 const std::vector<unsigned int> &ids = maps[c];
585 ocell->get_dof_indices(odofs);
588 cell_matrix =
typename Matrix::value_type();
590 for (
unsigned int i = 0;
591 i < space_dh.
get_fe().n_dofs_per_cell();
597 for (
unsigned int j = 0;
598 j < immersed_dh.
get_fe().n_dofs_per_cell();
601 const auto comp_j = immersed_dh.
get_fe()
604 if (space_gtl[comp_i] == immersed_gtl[comp_j])
605 for (
unsigned int oq = 0;
606 oq < o_fe_v.n_quadrature_points;
610 const unsigned int q = ids[oq];
613 (fe_v.shape_value(j, q) *
614 o_fe_v.shape_value(i, oq) * fe_v.JxW(q));
621 cell_matrix, odofs, immersed_constraints, dofs, matrix);
634 const double & epsilon,
648 ExcMessage(
"When epsilon is zero, you can only "
649 "call this function with dim1 <= dim0."));
664 const bool zero_is_distributed =
667 const bool one_is_distributed =
680 const bool outer_loop_on_zero =
681 (zero_is_distributed && !one_is_distributed) ||
685 const auto &fe0 = dh0.
get_fe();
686 const auto &fe1 = dh1.
get_fe();
689 std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
690 std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
692 if (outer_loop_on_zero)
698 std::vector<std::pair<
703 for (
const auto &cell0 :
706 intersection.resize(0);
710 boost::geometry::index::query(tree1,
711 boost::geometry::index::intersects(
713 std::back_inserter(intersection));
714 if (!intersection.empty())
716 cell0->get_dof_indices(dofs0);
717 for (
const auto &entry : intersection)
720 *entry.second, &dh1);
721 cell1->get_dof_indices(dofs1);
734 std::vector<std::pair<
739 for (
const auto &cell1 :
742 intersection.resize(0);
746 boost::geometry::index::query(tree0,
747 boost::geometry::index::intersects(
749 std::back_inserter(intersection));
750 if (!intersection.empty())
752 cell1->get_dof_indices(dofs1);
753 for (
const auto &entry : intersection)
756 *entry.second, &dh0);
757 cell0->get_dof_indices(dofs0);
769 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
773 const double & epsilon,
788 ExcMessage(
"When epsilon is zero, you can only "
789 "call this function with dim1 <= dim0."));
805 const bool zero_is_distributed =
808 const bool one_is_distributed =
821 const bool outer_loop_on_zero =
822 (zero_is_distributed && !one_is_distributed) ||
826 const auto &fe0 = dh0.
get_fe();
827 const auto &fe1 = dh1.
get_fe();
842 std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
843 std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
847 fe1.n_dofs_per_cell());
852 const auto >l0 = p.first;
853 const auto >l1 = p.second;
856 std::vector<double> kernel_values(quadrature1.
size());
858 auto assemble_one_pair = [&]() {
860 for (
unsigned int q0 = 0; q0 < quadrature0.
size(); ++q0)
863 kernel.
value_list(fev1.get_quadrature_points(), kernel_values);
864 for (
unsigned int j = 0; j < fe1.n_dofs_per_cell(); ++j)
866 const auto comp_j = fe1.system_to_component_index(j).first;
870 typename Matrix::value_type sum_q1 = {};
871 for (
unsigned int q1 = 0; q1 < quadrature1.
size(); ++q1)
873 fev1.shape_value(j, q1) * kernel_values[q1] * fev1.JxW(q1);
874 sum_q1 *= fev0.JxW(q0);
879 for (
unsigned int i = 0; i < fe0.n_dofs_per_cell(); ++i)
881 const auto comp_i = fe0.system_to_component_index(i).first;
883 gtl1[comp_j] == gtl0[comp_i])
884 cell_matrix(i, j) += fev0.shape_value(i, q0) * sum_q1;
895 if (outer_loop_on_zero)
901 std::vector<std::pair<
906 for (
const auto &cell0 :
909 intersection.resize(0);
913 boost::geometry::index::query(tree1,
914 boost::geometry::index::intersects(
916 std::back_inserter(intersection));
917 if (!intersection.empty())
919 cell0->get_dof_indices(dofs0);
921 for (
const auto &entry : intersection)
924 *entry.second, &dh1);
925 cell1->get_dof_indices(dofs1);
937 std::vector<std::pair<
942 for (
const auto &cell1 :
945 intersection.resize(0);
949 boost::geometry::index::query(tree0,
950 boost::geometry::index::intersects(
952 std::back_inserter(intersection));
953 if (!intersection.empty())
955 cell1->get_dof_indices(dofs1);
957 for (
const auto &entry : intersection)
960 *entry.second, &dh0);
961 cell0->get_dof_indices(dofs0);
970# include "coupling.inst"
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
void add_entries_local_to_global(const std::vector< size_type > &local_dof_indices, SparsityPatternType &sparsity_pattern, const bool keep_constrained_entries=true, const Table< 2, bool > &dof_mask=Table< 2, bool >()) const
void extend(const Number amount)
unsigned int size() const
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) const
types::global_dof_index n_dofs() const
unsigned int n_dofs_per_cell() const
unsigned int n_components() const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
virtual void value_list(const std::vector< Point< dim > > &points, std::vector< RangeNumberType > &values, const unsigned int component=0) const
virtual void set_radius(const double r)
virtual void set_center(const Point< dim > &p)
Abstract base class for mapping classes.
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
void reinit(const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > &cell)
unsigned int size() const
unsigned int n_active_cells() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_quadrature_points
Transformed quadrature points.
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
std::pair< std::vector< unsigned int >, std::vector< unsigned int > > compute_components_coupling(const ComponentMask &comps0, const ComponentMask &comps1, const FiniteElement< dim0, spacedim > &fe0, const FiniteElement< dim1, spacedim > &fe1)
std::pair< std::vector< Point< spacedim > >, std::vector< unsigned int > > qpoints_over_locally_owned_cells(const GridTools::Cache< dim0, spacedim > &cache, const DoFHandler< dim1, spacedim > &immersed_dh, const Quadrature< dim1 > &quad, const Mapping< dim1, spacedim > &immersed_mapping, const bool tria_is_parallel)
void create_coupling_sparsity_pattern(const DoFHandler< dim0, spacedim > &space_dh, const DoFHandler< dim1, spacedim > &immersed_dh, const Quadrature< dim1 > &quad, Sparsity &sparsity, const AffineConstraints< number > &constraints=AffineConstraints< number >(), const ComponentMask &space_comps=ComponentMask(), const ComponentMask &immersed_comps=ComponentMask(), const Mapping< dim0, spacedim > &space_mapping=StaticMappingQ1< dim0, spacedim >::mapping, const Mapping< dim1, spacedim > &immersed_mapping=StaticMappingQ1< dim1, spacedim >::mapping, const AffineConstraints< number > &immersed_constraints=AffineConstraints< number >())
void create_coupling_mass_matrix(const DoFHandler< dim0, spacedim > &space_dh, const DoFHandler< dim1, spacedim > &immersed_dh, const Quadrature< dim1 > &quad, Matrix &matrix, const AffineConstraints< typename Matrix::value_type > &constraints=AffineConstraints< typename Matrix::value_type >(), const ComponentMask &space_comps=ComponentMask(), const ComponentMask &immersed_comps=ComponentMask(), const Mapping< dim0, spacedim > &space_mapping=StaticMappingQ1< dim0, spacedim >::mapping, const Mapping< dim1, spacedim > &immersed_mapping=StaticMappingQ1< dim1, spacedim >::mapping, const AffineConstraints< typename Matrix::value_type > &immersed_constraints=AffineConstraints< typename Matrix::value_type >())
static const unsigned int invalid_unsigned_int