Reference documentation for deal.II version 9.3.2
\(\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\}}\)
coupling.cc
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2018 - 2021 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
17 #include <deal.II/base/point.h>
18 #include <deal.II/base/tensor.h>
19 
22 
23 #include <deal.II/fe/fe_values.h>
24 
27 
37 
39 
41 namespace NonMatching
42 {
43  namespace internal
44  {
62  template <int dim0, int dim1, int spacedim>
63  std::pair<std::vector<Point<spacedim>>, std::vector<unsigned int>>
66  const DoFHandler<dim1, spacedim> & immersed_dh,
67  const Quadrature<dim1> & quad,
68  const Mapping<dim1, spacedim> & immersed_mapping,
69  const bool tria_is_parallel)
70  {
71  const auto & immersed_fe = immersed_dh.get_fe();
72  std::vector<Point<spacedim>> points_over_local_cells;
73  // Keep track of which cells we actually used
74  std::vector<unsigned int> used_cells_ids;
75  {
76  FEValues<dim1, spacedim> fe_v(immersed_mapping,
77  immersed_fe,
78  quad,
80  unsigned int cell_id = 0;
81  for (const auto &cell : immersed_dh.active_cell_iterators())
82  {
83  bool use_cell = false;
84  if (tria_is_parallel)
85  {
86  const auto bbox = cell->bounding_box();
87  std::vector<std::pair<
90  out_vals;
91  cache.get_cell_bounding_boxes_rtree().query(
92  boost::geometry::index::intersects(bbox),
93  std::back_inserter(out_vals));
94  // Each bounding box corresponds to an active cell
95  // of the embedding triangulation: we now check if
96  // the current cell, of the embedded triangulation,
97  // overlaps a locally owned cell of the embedding one
98  for (const auto &bbox_it : out_vals)
99  if (bbox_it.second->is_locally_owned())
100  {
101  use_cell = true;
102  used_cells_ids.emplace_back(cell_id);
103  break;
104  }
105  }
106  else
107  // for sequential triangulations, simply use all cells
108  use_cell = true;
109 
110  if (use_cell)
111  {
112  // Reinitialize the cell and the fe_values
113  fe_v.reinit(cell);
114  const std::vector<Point<spacedim>> &x_points =
115  fe_v.get_quadrature_points();
116 
117  // Insert the points to the vector
118  points_over_local_cells.insert(points_over_local_cells.end(),
119  x_points.begin(),
120  x_points.end());
121  }
122  ++cell_id;
123  }
124  }
125  return {std::move(points_over_local_cells), std::move(used_cells_ids)};
126  }
127 
128 
135  template <int dim0, int dim1, int spacedim>
136  std::pair<std::vector<unsigned int>, std::vector<unsigned int>>
138  const ComponentMask & comps1,
141  {
142  // Take care of components
143  const ComponentMask mask0 =
144  (comps0.size() == 0 ? ComponentMask(fe0.n_components(), true) : comps0);
145 
146  const ComponentMask mask1 =
147  (comps1.size() == 0 ? ComponentMask(fe1.n_components(), true) : comps1);
148 
149  AssertDimension(mask0.size(), fe0.n_components());
150  AssertDimension(mask1.size(), fe1.n_components());
151 
152  // Global to local indices
153  std::vector<unsigned int> gtl0(fe0.n_components(),
155  std::vector<unsigned int> gtl1(fe1.n_components(),
157 
158  for (unsigned int i = 0, j = 0; i < gtl0.size(); ++i)
159  if (mask0[i])
160  gtl0[i] = j++;
161 
162  for (unsigned int i = 0, j = 0; i < gtl1.size(); ++i)
163  if (mask1[i])
164  gtl1[i] = j++;
165  return {gtl0, gtl1};
166  }
167  } // namespace internal
168 
169  template <int dim0,
170  int dim1,
171  int spacedim,
172  typename Sparsity,
173  typename number>
174  void
176  const DoFHandler<dim0, spacedim> &space_dh,
177  const DoFHandler<dim1, spacedim> &immersed_dh,
178  const Quadrature<dim1> & quad,
179  Sparsity & sparsity,
180  const AffineConstraints<number> & constraints,
181  const ComponentMask & space_comps,
182  const ComponentMask & immersed_comps,
183  const Mapping<dim0, spacedim> & space_mapping,
184  const Mapping<dim1, spacedim> & immersed_mapping)
185  {
187  space_mapping);
189  space_dh,
190  immersed_dh,
191  quad,
192  sparsity,
193  constraints,
194  space_comps,
195  immersed_comps,
196  immersed_mapping);
197  }
198 
199 
200 
201  template <int dim0,
202  int dim1,
203  int spacedim,
204  typename Sparsity,
205  typename number>
206  void
209  const DoFHandler<dim0, spacedim> & space_dh,
210  const DoFHandler<dim1, spacedim> & immersed_dh,
211  const Quadrature<dim1> & quad,
212  Sparsity & sparsity,
213  const AffineConstraints<number> & constraints,
214  const ComponentMask & space_comps,
215  const ComponentMask & immersed_comps,
216  const Mapping<dim1, spacedim> & immersed_mapping)
217  {
218  AssertDimension(sparsity.n_rows(), space_dh.n_dofs());
219  AssertDimension(sparsity.n_cols(), immersed_dh.n_dofs());
220  Assert(dim1 <= dim0,
221  ExcMessage("This function can only work if dim1 <= dim0"));
222  Assert((dynamic_cast<
224  &immersed_dh.get_triangulation()) == nullptr),
226 
227  const bool tria_is_parallel =
228  (dynamic_cast<const parallel::TriangulationBase<dim1, spacedim> *>(
229  &space_dh.get_triangulation()) != nullptr);
230  const auto &space_fe = space_dh.get_fe();
231  const auto &immersed_fe = immersed_dh.get_fe();
232 
233  // Dof indices
234  std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
235  std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
236 
237  // Take care of components
238  const ComponentMask space_c =
239  (space_comps.size() == 0 ? ComponentMask(space_fe.n_components(), true) :
240  space_comps);
241 
242  const ComponentMask immersed_c =
243  (immersed_comps.size() == 0 ?
244  ComponentMask(immersed_fe.n_components(), true) :
245  immersed_comps);
246 
247  AssertDimension(space_c.size(), space_fe.n_components());
248  AssertDimension(immersed_c.size(), immersed_fe.n_components());
249 
250  // Global to local indices
251  std::vector<unsigned int> space_gtl(space_fe.n_components(),
253  std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
255 
256  for (unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
257  if (space_c[i])
258  space_gtl[i] = j++;
259 
260  for (unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
261  if (immersed_c[i])
262  immersed_gtl[i] = j++;
263 
264  const unsigned int n_q_points = quad.size();
265  const unsigned int n_active_c =
266  immersed_dh.get_triangulation().n_active_cells();
267 
268  const auto qpoints_cells_data = internal::qpoints_over_locally_owned_cells(
269  cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
270 
271  const auto &points_over_local_cells = std::get<0>(qpoints_cells_data);
272  const auto &used_cells_ids = std::get<1>(qpoints_cells_data);
273 
274  // [TODO]: when the add_entries_local_to_global below will implement
275  // the version with the dof_mask, this should be uncommented.
276  //
277  // // Construct a dof_mask, used to distribute entries to the sparsity
278  // able< 2, bool > dof_mask(space_fe.n_dofs_per_cell(),
279  // immersed_fe.n_dofs_per_cell());
280  // of_mask.fill(false);
281  // or (unsigned int i=0; i<space_fe.n_dofs_per_cell(); ++i)
282  // {
283  // const auto comp_i = space_fe.system_to_component_index(i).first;
284  // if (space_gtl[comp_i] != numbers::invalid_unsigned_int)
285  // for (unsigned int j=0; j<immersed_fe.n_dofs_per_cell(); ++j)
286  // {
287  // const auto comp_j =
288  // immersed_fe.system_to_component_index(j).first; if
289  // (immersed_gtl[comp_j] == space_gtl[comp_i])
290  // dof_mask(i,j) = true;
291  // }
292  // }
293 
294 
295  // Get a list of outer cells, qpoints and maps.
296  const auto cpm =
297  GridTools::compute_point_locations(cache, points_over_local_cells);
298  const auto &all_cells = std::get<0>(cpm);
299  const auto &maps = std::get<2>(cpm);
300 
301  std::vector<
302  std::set<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
303  cell_sets(n_active_c);
304 
305  for (unsigned int i = 0; i < maps.size(); ++i)
306  {
307  // Quadrature points should be reasonably clustered:
308  // the following index keeps track of the last id
309  // where the current cell was inserted
310  unsigned int last_id = std::numeric_limits<unsigned int>::max();
311  unsigned int cell_id;
312  for (const unsigned int idx : maps[i])
313  {
314  // Find in which cell of immersed triangulation the point lies
315  if (tria_is_parallel)
316  cell_id = used_cells_ids[idx / n_q_points];
317  else
318  cell_id = idx / n_q_points;
319 
320  if (last_id != cell_id)
321  {
322  cell_sets[cell_id].insert(all_cells[i]);
323  last_id = cell_id;
324  }
325  }
326  }
327 
328  // Now we run on each cell of the immersed
329  // and build the sparsity
330  unsigned int i = 0;
331  for (const auto &cell : immersed_dh.active_cell_iterators())
332  {
333  // Reinitialize the cell
334  cell->get_dof_indices(dofs);
335 
336  // List of outer cells
337  const auto &cells = cell_sets[i];
338 
339  for (const auto &cell_c : cells)
340  {
341  // Get the ones in the current outer cell
342  typename DoFHandler<dim0, spacedim>::cell_iterator ocell(*cell_c,
343  &space_dh);
344  // Make sure we act only on locally_owned cells
345  if (ocell->is_locally_owned())
346  {
347  ocell->get_dof_indices(odofs);
348  // [TODO]: When the following function will be implemented
349  // for the case of non-trivial dof_mask, we should
350  // uncomment the missing part.
351  constraints.add_entries_local_to_global(
352  odofs, dofs, sparsity); //, true, dof_mask);
353  }
354  }
355  ++i;
356  }
357  }
358 
359 
360 
361  template <int dim0, int dim1, int spacedim, typename Matrix>
362  void
364  const DoFHandler<dim0, spacedim> & space_dh,
365  const DoFHandler<dim1, spacedim> & immersed_dh,
366  const Quadrature<dim1> & quad,
367  Matrix & matrix,
369  const ComponentMask & space_comps,
370  const ComponentMask & immersed_comps,
371  const Mapping<dim0, spacedim> & space_mapping,
372  const Mapping<dim1, spacedim> & immersed_mapping)
373  {
375  space_mapping);
377  space_dh,
378  immersed_dh,
379  quad,
380  matrix,
381  constraints,
382  space_comps,
383  immersed_comps,
384  immersed_mapping);
385  }
386 
387 
388 
389  template <int dim0, int dim1, int spacedim, typename Matrix>
390  void
392  const GridTools::Cache<dim0, spacedim> & cache,
393  const DoFHandler<dim0, spacedim> & space_dh,
394  const DoFHandler<dim1, spacedim> & immersed_dh,
395  const Quadrature<dim1> & quad,
396  Matrix & matrix,
398  const ComponentMask & space_comps,
399  const ComponentMask & immersed_comps,
400  const Mapping<dim1, spacedim> & immersed_mapping)
401  {
402  AssertDimension(matrix.m(), space_dh.n_dofs());
403  AssertDimension(matrix.n(), immersed_dh.n_dofs());
404  Assert(dim1 <= dim0,
405  ExcMessage("This function can only work if dim1 <= dim0"));
406  Assert((dynamic_cast<
408  &immersed_dh.get_triangulation()) == nullptr),
410 
411  const bool tria_is_parallel =
412  (dynamic_cast<const parallel::TriangulationBase<dim1, spacedim> *>(
413  &space_dh.get_triangulation()) != nullptr);
414 
415  const auto &space_fe = space_dh.get_fe();
416  const auto &immersed_fe = immersed_dh.get_fe();
417 
418  // Dof indices
419  std::vector<types::global_dof_index> dofs(immersed_fe.n_dofs_per_cell());
420  std::vector<types::global_dof_index> odofs(space_fe.n_dofs_per_cell());
421 
422  // Take care of components
423  const ComponentMask space_c =
424  (space_comps.size() == 0 ? ComponentMask(space_fe.n_components(), true) :
425  space_comps);
426 
427  const ComponentMask immersed_c =
428  (immersed_comps.size() == 0 ?
429  ComponentMask(immersed_fe.n_components(), true) :
430  immersed_comps);
431 
432  AssertDimension(space_c.size(), space_fe.n_components());
433  AssertDimension(immersed_c.size(), immersed_fe.n_components());
434 
435  std::vector<unsigned int> space_gtl(space_fe.n_components(),
437  std::vector<unsigned int> immersed_gtl(immersed_fe.n_components(),
439 
440  for (unsigned int i = 0, j = 0; i < space_gtl.size(); ++i)
441  if (space_c[i])
442  space_gtl[i] = j++;
443 
444  for (unsigned int i = 0, j = 0; i < immersed_gtl.size(); ++i)
445  if (immersed_c[i])
446  immersed_gtl[i] = j++;
447 
449  space_dh.get_fe().n_dofs_per_cell(),
450  immersed_dh.get_fe().n_dofs_per_cell());
451 
452  FEValues<dim1, spacedim> fe_v(immersed_mapping,
453  immersed_dh.get_fe(),
454  quad,
456  update_values);
457 
458  const unsigned int n_q_points = quad.size();
459  const unsigned int n_active_c =
460  immersed_dh.get_triangulation().n_active_cells();
461 
462  const auto used_cells_data = internal::qpoints_over_locally_owned_cells(
463  cache, immersed_dh, quad, immersed_mapping, tria_is_parallel);
464 
465  const auto &points_over_local_cells = std::get<0>(used_cells_data);
466  const auto &used_cells_ids = std::get<1>(used_cells_data);
467 
468  // Get a list of outer cells, qpoints and maps.
469  const auto cpm =
470  GridTools::compute_point_locations(cache, points_over_local_cells);
471  const auto &all_cells = std::get<0>(cpm);
472  const auto &all_qpoints = std::get<1>(cpm);
473  const auto &all_maps = std::get<2>(cpm);
474 
475  std::vector<
476  std::vector<typename Triangulation<dim0, spacedim>::active_cell_iterator>>
477  cell_container(n_active_c);
478  std::vector<std::vector<std::vector<Point<dim0>>>> qpoints_container(
479  n_active_c);
480  std::vector<std::vector<std::vector<unsigned int>>> maps_container(
481  n_active_c);
482 
483  // Cycle over all cells of underling mesh found
484  // call it omesh, elaborating the output
485  for (unsigned int o = 0; o < all_cells.size(); ++o)
486  {
487  for (unsigned int j = 0; j < all_maps[o].size(); ++j)
488  {
489  // Find the index of the "owner" cell and qpoint
490  // with regard to the immersed mesh
491  // Find in which cell of immersed triangulation the point lies
492  unsigned int cell_id;
493  if (tria_is_parallel)
494  cell_id = used_cells_ids[all_maps[o][j] / n_q_points];
495  else
496  cell_id = all_maps[o][j] / n_q_points;
497 
498  const unsigned int n_pt = all_maps[o][j] % n_q_points;
499 
500  // If there are no cells, we just add our data
501  if (cell_container[cell_id].empty())
502  {
503  cell_container[cell_id].emplace_back(all_cells[o]);
504  qpoints_container[cell_id].emplace_back(
505  std::vector<Point<dim0>>{all_qpoints[o][j]});
506  maps_container[cell_id].emplace_back(
507  std::vector<unsigned int>{n_pt});
508  }
509  // If there are already cells, we begin by looking
510  // at the last inserted cell, which is more likely:
511  else if (cell_container[cell_id].back() == all_cells[o])
512  {
513  qpoints_container[cell_id].back().emplace_back(
514  all_qpoints[o][j]);
515  maps_container[cell_id].back().emplace_back(n_pt);
516  }
517  else
518  {
519  // We don't need to check the last element
520  const auto cell_p = std::find(cell_container[cell_id].begin(),
521  cell_container[cell_id].end() - 1,
522  all_cells[o]);
523 
524  if (cell_p == cell_container[cell_id].end() - 1)
525  {
526  cell_container[cell_id].emplace_back(all_cells[o]);
527  qpoints_container[cell_id].emplace_back(
528  std::vector<Point<dim0>>{all_qpoints[o][j]});
529  maps_container[cell_id].emplace_back(
530  std::vector<unsigned int>{n_pt});
531  }
532  else
533  {
534  const unsigned int pos =
535  cell_p - cell_container[cell_id].begin();
536  qpoints_container[cell_id][pos].emplace_back(
537  all_qpoints[o][j]);
538  maps_container[cell_id][pos].emplace_back(n_pt);
539  }
540  }
541  }
542  }
543 
545  cell = immersed_dh.begin_active(),
546  endc = immersed_dh.end();
547 
548  for (unsigned int j = 0; cell != endc; ++cell, ++j)
549  {
550  // Reinitialize the cell and the fe_values
551  fe_v.reinit(cell);
552  cell->get_dof_indices(dofs);
553 
554  // Get a list of outer cells, qpoints and maps.
555  const auto &cells = cell_container[j];
556  const auto &qpoints = qpoints_container[j];
557  const auto &maps = maps_container[j];
558 
559  for (unsigned int c = 0; c < cells.size(); ++c)
560  {
561  // Get the ones in the current outer cell
563  *cells[c], &space_dh);
564  // Make sure we act only on locally_owned cells
565  if (ocell->is_locally_owned())
566  {
567  const std::vector<Point<dim0>> & qps = qpoints[c];
568  const std::vector<unsigned int> &ids = maps[c];
569 
570  FEValues<dim0, spacedim> o_fe_v(cache.get_mapping(),
571  space_dh.get_fe(),
572  qps,
573  update_values);
574  o_fe_v.reinit(ocell);
575  ocell->get_dof_indices(odofs);
576 
577  // Reset the matrices.
578  cell_matrix = typename Matrix::value_type();
579 
580  for (unsigned int i = 0;
581  i < space_dh.get_fe().n_dofs_per_cell();
582  ++i)
583  {
584  const auto comp_i =
585  space_dh.get_fe().system_to_component_index(i).first;
586  if (space_gtl[comp_i] != numbers::invalid_unsigned_int)
587  for (unsigned int j = 0;
588  j < immersed_dh.get_fe().n_dofs_per_cell();
589  ++j)
590  {
591  const auto comp_j = immersed_dh.get_fe()
593  .first;
594  if (space_gtl[comp_i] == immersed_gtl[comp_j])
595  for (unsigned int oq = 0;
596  oq < o_fe_v.n_quadrature_points;
597  ++oq)
598  {
599  // Get the corresponding q point
600  const unsigned int q = ids[oq];
601 
602  cell_matrix(i, j) +=
603  (fe_v.shape_value(j, q) *
604  o_fe_v.shape_value(i, oq) * fe_v.JxW(q));
605  }
606  }
607  }
608 
609  // Now assemble the matrices
611  odofs,
612  dofs,
613  matrix);
614  }
615  }
616  }
617  }
618 
619  template <int dim0,
620  int dim1,
621  int spacedim,
622  typename Sparsity,
623  typename Number>
624  void
626  const double & epsilon,
627  const GridTools::Cache<dim0, spacedim> &cache0,
628  const GridTools::Cache<dim1, spacedim> &cache1,
629  const DoFHandler<dim0, spacedim> & dh0,
630  const DoFHandler<dim1, spacedim> & dh1,
631  const Quadrature<dim1> & quad,
632  Sparsity & sparsity,
633  const AffineConstraints<Number> & constraints0,
634  const ComponentMask & comps0,
635  const ComponentMask & comps1)
636  {
637  if (epsilon == 0.0)
638  {
639  Assert(dim1 <= dim0,
640  ExcMessage("When epsilon is zero, you can only "
641  "call this function with dim1 <= dim0."));
643  dh0,
644  dh1,
645  quad,
646  sparsity,
647  constraints0,
648  comps0,
649  comps1,
650  cache1.get_mapping());
651  return;
652  }
653  AssertDimension(sparsity.n_rows(), dh0.n_dofs());
654  AssertDimension(sparsity.n_cols(), dh1.n_dofs());
655 
656  const bool zero_is_distributed =
658  *>(&dh0.get_triangulation()) != nullptr);
659  const bool one_is_distributed =
661  *>(&dh1.get_triangulation()) != nullptr);
662 
663  // We bail out if both are distributed triangulations
664  Assert(!zero_is_distributed || !one_is_distributed, ExcNotImplemented());
665 
666  // If we can loop on both, we decide where to make the outer loop according
667  // to the size of the triangulation. The reasoning is the following:
668  // - cost for accessing the tree: log(N)
669  // - cost for computing the intersection for each of the outer loop cells: M
670  // Total cost (besides the setup) is: M log(N)
671  // If we can, make sure M is the smaller number of the two.
672  const bool outer_loop_on_zero =
673  (zero_is_distributed && !one_is_distributed) ||
676 
677  const auto &fe0 = dh0.get_fe();
678  const auto &fe1 = dh1.get_fe();
679 
680  // Dof indices
681  std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
682  std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
683 
684  if (outer_loop_on_zero)
685  {
686  Assert(one_is_distributed == false, ExcInternalError());
687 
688  const auto &tree1 = cache1.get_cell_bounding_boxes_rtree();
689 
690  std::vector<std::pair<
693  intersection;
694 
695  for (const auto &cell0 : dh0.active_cell_iterators())
696  if (cell0->is_locally_owned())
697  {
698  intersection.resize(0);
699  BoundingBox<spacedim> box0 =
700  cache0.get_mapping().get_bounding_box(cell0);
701  box0.extend(epsilon);
702  boost::geometry::index::query(tree1,
703  boost::geometry::index::intersects(
704  box0),
705  std::back_inserter(intersection));
706  if (!intersection.empty())
707  {
708  cell0->get_dof_indices(dofs0);
709  for (const auto &entry : intersection)
710  {
712  *entry.second, &dh1);
713  cell1->get_dof_indices(dofs1);
714  constraints0.add_entries_local_to_global(dofs0,
715  dofs1,
716  sparsity);
717  }
718  }
719  }
720  }
721  else
722  {
723  Assert(zero_is_distributed == false, ExcInternalError());
724  const auto &tree0 = cache0.get_cell_bounding_boxes_rtree();
725 
726  std::vector<std::pair<
729  intersection;
730 
731  for (const auto &cell1 : dh1.active_cell_iterators())
732  if (cell1->is_locally_owned())
733  {
734  intersection.resize(0);
735  BoundingBox<spacedim> box1 =
736  cache1.get_mapping().get_bounding_box(cell1);
737  box1.extend(epsilon);
738  boost::geometry::index::query(tree0,
739  boost::geometry::index::intersects(
740  box1),
741  std::back_inserter(intersection));
742  if (!intersection.empty())
743  {
744  cell1->get_dof_indices(dofs1);
745  for (const auto &entry : intersection)
746  {
748  *entry.second, &dh0);
749  cell0->get_dof_indices(dofs0);
750  constraints0.add_entries_local_to_global(dofs0,
751  dofs1,
752  sparsity);
753  }
754  }
755  }
756  }
757  }
758 
759 
760 
761  template <int dim0, int dim1, int spacedim, typename Matrix>
762  void
765  const double & epsilon,
766  const GridTools::Cache<dim0, spacedim> & cache0,
767  const GridTools::Cache<dim1, spacedim> & cache1,
768  const DoFHandler<dim0, spacedim> & dh0,
769  const DoFHandler<dim1, spacedim> & dh1,
770  const Quadrature<dim0> & quadrature0,
771  const Quadrature<dim1> & quadrature1,
772  Matrix & matrix,
774  const ComponentMask & comps0,
775  const ComponentMask & comps1)
776  {
777  if (epsilon == 0)
778  {
779  Assert(dim1 <= dim0,
780  ExcMessage("When epsilon is zero, you can only "
781  "call this function with dim1 <= dim0."));
783  dh0,
784  dh1,
785  quadrature1,
786  matrix,
787  constraints0,
788  comps0,
789  comps1,
790  cache1.get_mapping());
791  return;
792  }
793 
794  AssertDimension(matrix.m(), dh0.n_dofs());
795  AssertDimension(matrix.n(), dh1.n_dofs());
796 
797  const bool zero_is_distributed =
799  *>(&dh0.get_triangulation()) != nullptr);
800  const bool one_is_distributed =
802  *>(&dh1.get_triangulation()) != nullptr);
803 
804  // We bail out if both are distributed triangulations
805  Assert(!zero_is_distributed || !one_is_distributed, ExcNotImplemented());
806 
807  // If we can loop on both, we decide where to make the outer loop according
808  // to the size of the triangulation. The reasoning is the following:
809  // - cost for accessing the tree: log(N)
810  // - cost for computing the intersection for each of the outer loop cells: M
811  // Total cost (besides the setup) is: M log(N)
812  // If we can, make sure M is the smaller number of the two.
813  const bool outer_loop_on_zero =
814  (zero_is_distributed && !one_is_distributed) ||
817 
818  const auto &fe0 = dh0.get_fe();
819  const auto &fe1 = dh1.get_fe();
820 
822  fe0,
823  quadrature0,
826 
828  fe1,
829  quadrature1,
832 
833  // Dof indices
834  std::vector<types::global_dof_index> dofs0(fe0.n_dofs_per_cell());
835  std::vector<types::global_dof_index> dofs1(fe1.n_dofs_per_cell());
836 
837  // Local Matrix
839  fe1.n_dofs_per_cell());
840 
841  // Global to local indices
842  const auto p =
843  internal::compute_components_coupling(comps0, comps1, fe0, fe1);
844  const auto &gtl0 = p.first;
845  const auto &gtl1 = p.second;
846 
847  kernel.set_radius(epsilon);
848  std::vector<double> kernel_values(quadrature1.size());
849 
850  auto assemble_one_pair = [&]() {
851  cell_matrix = 0;
852  for (unsigned int q0 = 0; q0 < quadrature0.size(); ++q0)
853  {
854  kernel.set_center(fev0.quadrature_point(q0));
855  kernel.value_list(fev1.get_quadrature_points(), kernel_values);
856  for (unsigned int j = 0; j < fe1.n_dofs_per_cell(); ++j)
857  {
858  const auto comp_j = fe1.system_to_component_index(j).first;
859 
860  // First compute the part of the integral that does not
861  // depend on i
862  typename Matrix::value_type sum_q1 = {};
863  for (unsigned int q1 = 0; q1 < quadrature1.size(); ++q1)
864  sum_q1 +=
865  fev1.shape_value(j, q1) * kernel_values[q1] * fev1.JxW(q1);
866  sum_q1 *= fev0.JxW(q0);
867 
868  // Now compute the main integral with the sum over q1 already
869  // completed - this gives a cubic complexity as usual rather
870  // than a quartic one with naive loops
871  for (unsigned int i = 0; i < fe0.n_dofs_per_cell(); ++i)
872  {
873  const auto comp_i = fe0.system_to_component_index(i).first;
874  if (gtl0[comp_i] != numbers::invalid_unsigned_int &&
875  gtl1[comp_j] == gtl0[comp_i])
876  cell_matrix(i, j) += fev0.shape_value(i, q0) * sum_q1;
877  }
878  }
879  }
880 
882  dofs0,
883  dofs1,
884  matrix);
885  };
886 
887  if (outer_loop_on_zero)
888  {
889  Assert(one_is_distributed == false, ExcInternalError());
890 
891  const auto &tree1 = cache1.get_cell_bounding_boxes_rtree();
892 
893  std::vector<std::pair<
896  intersection;
897 
898  for (const auto &cell0 : dh0.active_cell_iterators())
899  if (cell0->is_locally_owned())
900  {
901  intersection.resize(0);
902  BoundingBox<spacedim> box0 =
903  cache0.get_mapping().get_bounding_box(cell0);
904  box0.extend(epsilon);
905  boost::geometry::index::query(tree1,
906  boost::geometry::index::intersects(
907  box0),
908  std::back_inserter(intersection));
909  if (!intersection.empty())
910  {
911  cell0->get_dof_indices(dofs0);
912  fev0.reinit(cell0);
913  for (const auto &entry : intersection)
914  {
916  *entry.second, &dh1);
917  cell1->get_dof_indices(dofs1);
918  fev1.reinit(cell1);
919  assemble_one_pair();
920  }
921  }
922  }
923  }
924  else
925  {
926  Assert(zero_is_distributed == false, ExcInternalError());
927  const auto &tree0 = cache0.get_cell_bounding_boxes_rtree();
928 
929  std::vector<std::pair<
932  intersection;
933 
934  for (const auto &cell1 : dh1.active_cell_iterators())
935  if (cell1->is_locally_owned())
936  {
937  intersection.resize(0);
938  BoundingBox<spacedim> box1 =
939  cache1.get_mapping().get_bounding_box(cell1);
940  box1.extend(epsilon);
941  boost::geometry::index::query(tree0,
942  boost::geometry::index::intersects(
943  box1),
944  std::back_inserter(intersection));
945  if (!intersection.empty())
946  {
947  cell1->get_dof_indices(dofs1);
948  fev1.reinit(cell1);
949  for (const auto &entry : intersection)
950  {
952  *entry.second, &dh0);
953  cell0->get_dof_indices(dofs0);
954  fev0.reinit(cell0);
955  assemble_one_pair();
956  }
957  }
958  }
959  }
960  }
961 
962 #include "coupling.inst"
963 } // namespace NonMatching
964 
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 Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) const
types::global_dof_index n_dofs() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const double & shape_value(const unsigned int function_no, const unsigned int point_no) const
const unsigned int n_quadrature_points
Definition: fe_values.h:2432
const Point< spacedim > & quadrature_point(const unsigned int q) const
double JxW(const unsigned int quadrature_point) const
const std::vector< Point< spacedim > > & get_quadrature_points() const
void reinit(const TriaIterator< DoFCellAccessor< dim, spacedim, level_dof_access >> &cell)
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)
const Mapping< dim, spacedim > & get_mapping() const
const RTree< std::pair< BoundingBox< spacedim >, typename Triangulation< dim, spacedim >::active_cell_iterator > > & get_cell_bounding_boxes_rtree() const
Abstract base class for mapping classes.
Definition: mapping.h:304
virtual BoundingBox< spacedim > get_bounding_box(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition: point.h:111
unsigned int size() const
unsigned int n_active_cells() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
@ 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 & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcNotImplemented()
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
static ::ExceptionBase & ExcMessage(std::string arg1)
typename ActiveSelector::cell_iterator cell_iterator
Definition: dof_handler.h:466
typename ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:438
return_type compute_point_locations(const Cache< dim, spacedim > &cache, const std::vector< Point< spacedim >> &points, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell_hint=typename Triangulation< dim, spacedim >::active_cell_iterator())
Definition: grid_tools.cc:5425
@ matrix
Contents is actually a matrix.
void cell_matrix(FullMatrix< double > &M, const FEValuesBase< dim > &fe, const FEValuesBase< dim > &fetest, const ArrayView< const std::vector< double >> &velocity, const double factor=1.)
Definition: advection.h:75
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)
Definition: coupling.cc:137
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)
Definition: coupling.cc:64
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)
Definition: coupling.cc:175
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)
Definition: coupling.cc:363
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
VectorType::value_type * begin(VectorType &V)
VectorType::value_type * end(VectorType &V)
static const unsigned int invalid_unsigned_int
Definition: types.h:196