64 template <
int dim0,
int dim1,
int spacedim>
65 std::pair<std::vector<Point<spacedim>>, std::vector<unsigned int>>
71 const bool tria_is_parallel)
73 const auto & immersed_fe = immersed_dh.
get_fe();
74 std::vector<Point<spacedim>> points_over_local_cells;
76 std::vector<unsigned int> used_cells_ids;
82 unsigned int cell_id = 0;
85 bool use_cell =
false;
88 const auto bbox = cell->bounding_box();
89 std::vector<std::pair<
94 boost::geometry::index::intersects(bbox),
95 std::back_inserter(out_vals));
100 for (
const auto &bbox_it : out_vals)
101 if (bbox_it.second->is_locally_owned())
104 used_cells_ids.emplace_back(cell_id);
116 const std::vector<Point<spacedim>> &x_points =
117 fe_v.get_quadrature_points();
120 points_over_local_cells.insert(points_over_local_cells.end(),
127 return {std::move(points_over_local_cells), std::move(used_cells_ids)};
137 template <
int dim0,
int dim1,
int spacedim>
138 std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
160 for (
unsigned int i = 0, j = 0; i < gtl0.size(); ++i)
164 for (
unsigned int i = 0, j = 0; i < gtl1.size(); ++i)
171 template <
int dim0,
int dim1,
int spacedim,
typename number>
196 immersed_constraints);
201 template <
int dim0,
int dim1,
int spacedim,
typename number>
218 ExcMessage(
"This function can only work if dim1 <= dim0"));
224 const bool tria_is_parallel =
227 const auto &space_fe = space_dh.
get_fe();
228 const auto &immersed_fe = immersed_dh.
get_fe();
231 std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
232 std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
240 (immersed_comps.
size() == 0 ?
248 std::vector<unsigned int> space_gtl(space_fe.n_components(),
250 std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
253 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
257 for (
unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
259 immersed_gtl[i] = j++;
261 const unsigned int n_q_points = quad.
size();
262 const unsigned int n_active_c =
266 cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
268 const auto &points_over_local_cells = std::get<0>(qpoints_cells_data);
269 const auto &used_cells_ids = std::get<1>(qpoints_cells_data);
295 const auto &all_cells = std::get<0>(cpm);
296 const auto &maps = std::get<2>(cpm);
299 std::set<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
300 cell_sets(n_active_c);
302 for (
unsigned int i = 0; i < maps.size(); ++i)
307 unsigned int last_id = std::numeric_limits<unsigned int>::max();
308 unsigned int cell_id;
309 for (
const unsigned int idx : maps[i])
312 if (tria_is_parallel)
313 cell_id = used_cells_ids[idx / n_q_points];
315 cell_id = idx / n_q_points;
317 if (last_id != cell_id)
319 cell_sets[cell_id].insert(all_cells[i]);
331 cell->get_dof_indices(dofs);
334 const auto &cells = cell_sets[i];
336 for (
const auto &cell_c : cells)
342 if (ocell->is_locally_owned())
344 ocell->get_dof_indices(odofs);
350 immersed_constraints,
361 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
386 immersed_constraints);
391 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
408 ExcMessage(
"This function can only work if dim1 <= dim0"));
414 const bool tria_is_parallel =
418 const auto &space_fe = space_dh.
get_fe();
419 const auto &immersed_fe = immersed_dh.
get_fe();
422 std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
423 std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
431 (immersed_comps.
size() == 0 ?
438 std::vector<unsigned int> space_gtl(space_fe.n_components(),
440 std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
443 for (
unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
447 for (
unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
449 immersed_gtl[i] = j++;
461 const unsigned int n_q_points = quad.
size();
462 const unsigned int n_active_c =
466 cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
468 const auto &points_over_local_cells = std::get<0>(used_cells_data);
469 const auto &used_cells_ids = std::get<1>(used_cells_data);
474 const auto &all_cells = std::get<0>(cpm);
475 const auto &all_qpoints = std::get<1>(cpm);
476 const auto &all_maps = std::get<2>(cpm);
479 std::vector<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
480 cell_container(n_active_c);
481 std::vector<std::vector<std::vector<Point<dim0>>>> qpoints_container(
483 std::vector<std::vector<std::vector<unsigned int>>> maps_container(
488 for (
unsigned int o = 0; o < all_cells.size(); ++o)
490 for (
unsigned int j = 0; j < all_maps[o].size(); ++j)
495 unsigned int cell_id;
496 if (tria_is_parallel)
497 cell_id = used_cells_ids[all_maps[o][j] / n_q_points];
499 cell_id = all_maps[o][j] / n_q_points;
501 const unsigned int n_pt = all_maps[o][j] % n_q_points;
504 if (cell_container[cell_id].empty())
506 cell_container[cell_id].emplace_back(all_cells[o]);
507 qpoints_container[cell_id].emplace_back(
509 maps_container[cell_id].emplace_back(
510 std::vector<unsigned int>{n_pt});
514 else if (cell_container[cell_id].back() == all_cells[o])
516 qpoints_container[cell_id].back().emplace_back(
518 maps_container[cell_id].back().emplace_back(n_pt);
523 const auto cell_p = std::find(cell_container[cell_id].begin(),
524 cell_container[cell_id].end() - 1,
527 if (cell_p == cell_container[cell_id].end() - 1)
529 cell_container[cell_id].emplace_back(all_cells[o]);
530 qpoints_container[cell_id].emplace_back(
532 maps_container[cell_id].emplace_back(
533 std::vector<unsigned int>{n_pt});
537 const unsigned int pos =
538 cell_p - cell_container[cell_id].begin();
539 qpoints_container[cell_id][pos].emplace_back(
541 maps_container[cell_id][pos].emplace_back(n_pt);
549 endc = immersed_dh.
end();
551 for (
unsigned int j = 0; cell != endc; ++cell, ++j)
555 cell->get_dof_indices(dofs);
558 const auto &cells = cell_container[j];
559 const auto &qpoints = qpoints_container[j];
560 const auto &maps = maps_container[j];
562 for (
unsigned int c = 0; c < cells.size(); ++c)
566 *cells[c], &space_dh);
568 if (ocell->is_locally_owned())
570 const std::vector<Point<dim0>> & qps = qpoints[c];
571 const std::vector<unsigned int> &ids = maps[c];
578 ocell->get_dof_indices(odofs);
581 cell_matrix =
typename Matrix::value_type();
583 for (
unsigned int i = 0;
584 i < space_dh.
get_fe().n_dofs_per_cell();
590 for (
unsigned int j = 0;
591 j < immersed_dh.
get_fe().n_dofs_per_cell();
594 const auto comp_j = immersed_dh.
get_fe()
597 if (space_gtl[comp_i] == immersed_gtl[comp_j])
598 for (
unsigned int oq = 0;
599 oq < o_fe_v.n_quadrature_points;
603 const unsigned int q = ids[oq];
606 (fe_v.shape_value(j, q) *
607 o_fe_v.shape_value(i, oq) * fe_v.JxW(q));
614 cell_matrix, odofs, immersed_constraints, dofs, matrix);
620 template <
int dim0,
int dim1,
int spacedim,
typename Number>
623 const double & epsilon,
637 ExcMessage(
"When epsilon is zero, you can only "
638 "call this function with dim1 <= dim0."));
653 const bool zero_is_distributed =
656 const bool one_is_distributed =
669 const bool outer_loop_on_zero =
670 (zero_is_distributed && !one_is_distributed) ||
674 const auto &fe0 = dh0.
get_fe();
675 const auto &fe1 = dh1.
get_fe();
678 std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
679 std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
681 if (outer_loop_on_zero)
687 std::vector<std::pair<
692 for (
const auto &cell0 :
695 intersection.resize(0);
699 boost::geometry::index::query(tree1,
700 boost::geometry::index::intersects(
702 std::back_inserter(intersection));
703 if (!intersection.empty())
705 cell0->get_dof_indices(dofs0);
706 for (
const auto &entry : intersection)
709 *entry.second, &dh1);
710 cell1->get_dof_indices(dofs1);
723 std::vector<std::pair<
728 for (
const auto &cell1 :
731 intersection.resize(0);
735 boost::geometry::index::query(tree0,
736 boost::geometry::index::intersects(
738 std::back_inserter(intersection));
739 if (!intersection.empty())
741 cell1->get_dof_indices(dofs1);
742 for (
const auto &entry : intersection)
745 *entry.second, &dh0);
746 cell0->get_dof_indices(dofs0);
758 template <
int dim0,
int dim1,
int spacedim,
typename Matrix>
762 const double & epsilon,
777 ExcMessage(
"When epsilon is zero, you can only "
778 "call this function with dim1 <= dim0."));
794 const bool zero_is_distributed =
797 const bool one_is_distributed =
810 const bool outer_loop_on_zero =
811 (zero_is_distributed && !one_is_distributed) ||
815 const auto &fe0 = dh0.
get_fe();
816 const auto &fe1 = dh1.
get_fe();
831 std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
832 std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
836 fe1.n_dofs_per_cell());
841 const auto >l0 = p.first;
842 const auto >l1 = p.second;
845 std::vector<double> kernel_values(quadrature1.
size());
847 auto assemble_one_pair = [&]() {
849 for (
unsigned int q0 = 0; q0 < quadrature0.
size(); ++q0)
852 kernel.
value_list(fev1.get_quadrature_points(), kernel_values);
853 for (
unsigned int j = 0; j < fe1.n_dofs_per_cell(); ++j)
855 const auto comp_j = fe1.system_to_component_index(j).first;
859 typename Matrix::value_type sum_q1 = {};
860 for (
unsigned int q1 = 0; q1 < quadrature1.
size(); ++q1)
862 fev1.shape_value(j, q1) * kernel_values[q1] * fev1.JxW(q1);
863 sum_q1 *= fev0.JxW(q0);
868 for (
unsigned int i = 0; i < fe0.n_dofs_per_cell(); ++i)
870 const auto comp_i = fe0.system_to_component_index(i).first;
872 gtl1[comp_j] == gtl0[comp_i])
873 cell_matrix(i, j) += fev0.shape_value(i, q0) * sum_q1;
884 if (outer_loop_on_zero)
890 std::vector<std::pair<
895 for (
const auto &cell0 :
898 intersection.resize(0);
902 boost::geometry::index::query(tree1,
903 boost::geometry::index::intersects(
905 std::back_inserter(intersection));
906 if (!intersection.empty())
908 cell0->get_dof_indices(dofs0);
910 for (
const auto &entry : intersection)
913 *entry.second, &dh1);
914 cell1->get_dof_indices(dofs1);
926 std::vector<std::pair<
931 for (
const auto &cell1 :
934 intersection.resize(0);
938 boost::geometry::index::query(tree0,
939 boost::geometry::index::intersects(
941 std::back_inserter(intersection));
942 if (!intersection.empty())
944 cell1->get_dof_indices(dofs1);
946 for (
const auto &entry : intersection)
949 *entry.second, &dh0);
950 cell0->get_dof_indices(dofs0);
959# 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)
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_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 >())
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=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 >())
static const unsigned int invalid_unsigned_int