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