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)
170 template <
int dim0,
int dim1,
int spacedim,
typename number>
195 immersed_constraints);
200 template <
int dim0,
int dim1,
int spacedim,
typename number>
217 ExcMessage(
"This function can only work if dim1 <= dim0"));
223 const bool tria_is_parallel =
226 const auto &space_fe = space_dh.
get_fe();
227 const auto &immersed_fe = immersed_dh.
get_fe();
230 std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
231 std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
239 (immersed_comps.
size() == 0 ?
247 std::vector<unsigned int> space_gtl(space_fe.n_components(),
249 std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
252 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
256 for (
unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
258 immersed_gtl[i] = j++;
260 const unsigned int n_q_points = quad.
size();
261 const unsigned int n_active_c =
265 cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
267 const auto &points_over_local_cells = std::get<0>(qpoints_cells_data);
268 const auto &used_cells_ids = std::get<1>(qpoints_cells_data);
294 const auto &all_cells = std::get<0>(cpm);
295 const auto &maps = std::get<2>(cpm);
298 std::set<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
299 cell_sets(n_active_c);
301 for (
unsigned int i = 0; i < maps.size(); ++i)
306 unsigned int last_id = std::numeric_limits<unsigned int>::max();
307 unsigned int cell_id;
308 for (
const unsigned int idx : maps[i])
311 if (tria_is_parallel)
312 cell_id = used_cells_ids[idx / n_q_points];
314 cell_id = idx / n_q_points;
316 if (last_id != cell_id)
318 cell_sets[cell_id].insert(all_cells[i]);
330 cell->get_dof_indices(dofs);
333 const auto &cells = cell_sets[i];
335 for (
const auto &cell_c : cells)
341 if (ocell->is_locally_owned())
343 ocell->get_dof_indices(odofs);
349 immersed_constraints,
360 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
385 immersed_constraints);
390 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
407 ExcMessage(
"This function can only work if dim1 <= dim0"));
413 const bool tria_is_parallel =
417 const auto &space_fe = space_dh.
get_fe();
418 const auto &immersed_fe = immersed_dh.
get_fe();
421 std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
422 std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
430 (immersed_comps.
size() == 0 ?
437 std::vector<unsigned int> space_gtl(space_fe.n_components(),
439 std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
442 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
446 for (
unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
448 immersed_gtl[i] = j++;
460 const unsigned int n_q_points = quad.
size();
461 const unsigned int n_active_c =
465 cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
467 const auto &points_over_local_cells = std::get<0>(used_cells_data);
468 const auto &used_cells_ids = std::get<1>(used_cells_data);
473 const auto &all_cells = std::get<0>(cpm);
474 const auto &all_qpoints = std::get<1>(cpm);
475 const auto &all_maps = std::get<2>(cpm);
478 std::vector<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
479 cell_container(n_active_c);
480 std::vector<std::vector<std::vector<Point<dim0>>>> qpoints_container(
482 std::vector<std::vector<std::vector<unsigned int>>> maps_container(
487 for (
unsigned int o = 0; o < all_cells.size(); ++o)
489 for (
unsigned int j = 0; j < all_maps[o].size(); ++j)
494 unsigned int cell_id;
495 if (tria_is_parallel)
496 cell_id = used_cells_ids[all_maps[o][j] / n_q_points];
498 cell_id = all_maps[o][j] / n_q_points;
500 const unsigned int n_pt = all_maps[o][j] % n_q_points;
503 if (cell_container[cell_id].empty())
505 cell_container[cell_id].emplace_back(all_cells[o]);
506 qpoints_container[cell_id].emplace_back(
508 maps_container[cell_id].emplace_back(
509 std::vector<unsigned int>{n_pt});
513 else if (cell_container[cell_id].back() == all_cells[o])
515 qpoints_container[cell_id].back().emplace_back(
517 maps_container[cell_id].back().emplace_back(n_pt);
522 const auto cell_p = std::find(cell_container[cell_id].begin(),
523 cell_container[cell_id].end() - 1,
526 if (cell_p == cell_container[cell_id].end() - 1)
528 cell_container[cell_id].emplace_back(all_cells[o]);
529 qpoints_container[cell_id].emplace_back(
531 maps_container[cell_id].emplace_back(
532 std::vector<unsigned int>{n_pt});
536 const unsigned int pos =
537 cell_p - cell_container[cell_id].begin();
538 qpoints_container[cell_id][pos].emplace_back(
540 maps_container[cell_id][pos].emplace_back(n_pt);
548 endc = immersed_dh.
end();
550 for (
unsigned int j = 0; cell != endc; ++cell, ++j)
554 cell->get_dof_indices(dofs);
557 const auto &cells = cell_container[j];
558 const auto &qpoints = qpoints_container[j];
559 const auto &maps = maps_container[j];
561 for (
unsigned int c = 0; c < cells.size(); ++c)
565 *cells[c], &space_dh);
567 if (ocell->is_locally_owned())
569 const std::vector<Point<dim0>> &qps = qpoints[c];
570 const std::vector<unsigned int> &ids = maps[c];
577 ocell->get_dof_indices(odofs);
580 cell_matrix =
typename Matrix::value_type();
582 for (
unsigned int i = 0;
583 i < space_dh.
get_fe().n_dofs_per_cell();
589 for (
unsigned int j = 0;
590 j < immersed_dh.
get_fe().n_dofs_per_cell();
593 const auto comp_j = immersed_dh.
get_fe()
596 if (space_gtl[comp_i] == immersed_gtl[comp_j])
597 for (
unsigned int oq = 0;
598 oq < o_fe_v.n_quadrature_points;
602 const unsigned int q = ids[oq];
605 (fe_v.shape_value(j, q) *
606 o_fe_v.shape_value(i, oq) * fe_v.JxW(q));
613 cell_matrix, odofs, immersed_constraints, dofs, matrix);
619 template <
int dim0,
int dim1,
int spacedim,
typename Number>
622 const double &epsilon,
636 ExcMessage(
"When epsilon is zero, you can only "
637 "call this function with dim1 <= dim0."));
652 const bool zero_is_distributed =
655 const bool one_is_distributed =
668 const bool outer_loop_on_zero =
669 (zero_is_distributed && !one_is_distributed) ||
673 const auto &fe0 = dh0.
get_fe();
674 const auto &fe1 = dh1.
get_fe();
677 std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
678 std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
680 if (outer_loop_on_zero)
686 std::vector<std::pair<
691 for (
const auto &cell0 :
694 intersection.resize(0);
698 boost::geometry::index::query(tree1,
699 boost::geometry::index::intersects(
701 std::back_inserter(intersection));
702 if (!intersection.empty())
704 cell0->get_dof_indices(dofs0);
705 for (
const auto &entry : intersection)
708 *entry.second, &dh1);
709 cell1->get_dof_indices(dofs1);
722 std::vector<std::pair<
727 for (
const auto &cell1 :
730 intersection.resize(0);
734 boost::geometry::index::query(tree0,
735 boost::geometry::index::intersects(
737 std::back_inserter(intersection));
738 if (!intersection.empty())
740 cell1->get_dof_indices(dofs1);
741 for (
const auto &entry : intersection)
744 *entry.second, &dh0);
745 cell0->get_dof_indices(dofs0);
757 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
761 const double &epsilon,
776 ExcMessage(
"When epsilon is zero, you can only "
777 "call this function with dim1 <= dim0."));
793 const bool zero_is_distributed =
796 const bool one_is_distributed =
809 const bool outer_loop_on_zero =
810 (zero_is_distributed && !one_is_distributed) ||
814 const auto &fe0 = dh0.
get_fe();
815 const auto &fe1 = dh1.
get_fe();
830 std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
831 std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
835 fe1.n_dofs_per_cell());
840 const auto >l0 = p.first;
841 const auto >l1 = p.second;
844 std::vector<double> kernel_values(quadrature1.
size());
846 auto assemble_one_pair = [&]() {
848 for (
unsigned int q0 = 0; q0 < quadrature0.
size(); ++q0)
851 kernel.
value_list(fev1.get_quadrature_points(), kernel_values);
852 for (
unsigned int j = 0; j < fe1.n_dofs_per_cell(); ++j)
854 const auto comp_j = fe1.system_to_component_index(j).first;
858 typename Matrix::value_type sum_q1 = {};
859 for (
unsigned int q1 = 0; q1 < quadrature1.
size(); ++q1)
861 fev1.shape_value(j, q1) * kernel_values[q1] * fev1.JxW(q1);
862 sum_q1 *= fev0.JxW(q0);
867 for (
unsigned int i = 0; i < fe0.n_dofs_per_cell(); ++i)
869 const auto comp_i = fe0.system_to_component_index(i).first;
871 gtl1[comp_j] == gtl0[comp_i])
872 cell_matrix(i, j) += fev0.shape_value(i, q0) * sum_q1;
883 if (outer_loop_on_zero)
889 std::vector<std::pair<
894 for (
const auto &cell0 :
897 intersection.resize(0);
901 boost::geometry::index::query(tree1,
902 boost::geometry::index::intersects(
904 std::back_inserter(intersection));
905 if (!intersection.empty())
907 cell0->get_dof_indices(dofs0);
909 for (
const auto &entry : intersection)
912 *entry.second, &dh1);
913 cell1->get_dof_indices(dofs1);
925 std::vector<std::pair<
930 for (
const auto &cell1 :
933 intersection.resize(0);
937 boost::geometry::index::query(tree0,
938 boost::geometry::index::intersects(
940 std::back_inserter(intersection));
941 if (!intersection.empty())
943 cell1->get_dof_indices(dofs1);
945 for (
const auto &entry : intersection)
948 *entry.second, &dh0);
949 cell0->get_dof_indices(dofs0);
958# 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, SparsityPatternBase &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 types::fe_index 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, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int)
unsigned int size() const
unsigned int n_active_cells() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
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
@ update_values
Shape function values.
@ update_JxW_values
Transformed quadrature weights.
@ update_quadrature_points
Transformed quadrature points.
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, SparsityPatternBase &sparsity, const AffineConstraints< number > &constraints={}, const ComponentMask &space_comps={}, const ComponentMask &immersed_comps={}, 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={}, const ComponentMask &immersed_comps={}, 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