PolyDEAL
 
Loading...
Searching...
No Matches
agglomeration_handler.cc
Go to the documentation of this file.
1// -----------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception OR LGPL-2.1-or-later
4// Copyright (C) XXXX - YYYY by the polyDEAL authors
5//
6// This file is part of the polyDEAL library.
7//
8// Detailed license information governing the source code
9// can be found in LICENSE.md at the top level directory.
10//
11// -----------------------------------------------------------------------------
12
13
14
16
18#include <fe_agglodgp.h>
19
20template <int dim, int spacedim>
22 const GridTools::Cache<dim, spacedim> &cache_tria)
23 : cached_tria(std::make_unique<GridTools::Cache<dim, spacedim>>(
24 cache_tria.get_triangulation(),
25 cache_tria.get_mapping()))
26 , communicator(cache_tria.get_triangulation().get_mpi_communicator())
27{
28 Assert(dim == spacedim, ExcNotImplemented("Not available with codim > 0"));
29 Assert(dim == 2 || dim == 3, ExcImpossibleInDim(1));
31 &cached_tria->get_triangulation()) == nullptr),
32 ExcNotImplemented());
33 Assert(cached_tria->get_triangulation().n_active_cells() > 0,
34 ExcMessage(
35 "The triangulation must not be empty upon calling this function."));
36
38 hybrid_mesh = false;
40}
41
42
43
44template <int dim, int spacedim>
47 const AgglomerationContainer &cells)
48{
49 Assert(cells.size() > 0, ExcMessage("No cells to be agglomerated."));
50
51 if (cells.size() == 1)
52 hybrid_mesh = true; // mesh is made also by classical cells
53
54 // First index drives the selection of the master cell. After that, store the
55 // master cell.
56 const types::global_cell_index global_master_idx =
57 cells[0]->global_active_cell_index();
58 const types::global_cell_index master_idx = cells[0]->active_cell_index();
59 master_cells_container.push_back(cells[0]);
60 master_slave_relationships[global_master_idx] = -1;
61
62 const typename DoFHandler<dim>::active_cell_iterator cell_dh =
63 cells[0]->as_dof_handler_iterator(agglo_dh);
64 cell_dh->set_active_fe_index(CellAgglomerationType::master);
65
66 // Store slave cells and save the relationship with the parent
67 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
68 slaves;
69 slaves.reserve(cells.size() - 1);
70 // exclude first cell since it's the master cell
71 for (auto it = ++cells.begin(); it != cells.end(); ++it)
72 {
73 slaves.push_back(*it);
74 master_slave_relationships[(*it)->global_active_cell_index()] =
75 global_master_idx; // mark each slave
76 master_slave_relationships_iterators[(*it)->active_cell_index()] =
77 cells[0];
78
79 const typename DoFHandler<dim>::active_cell_iterator cell =
80 (*it)->as_dof_handler_iterator(agglo_dh);
81 cell->set_active_fe_index(CellAgglomerationType::slave); // slave cell
82
83 // If we have a p::d::T, check that all cells are in the same subdomain.
84 // If serial, just check that the subdomain_id is invalid.
85 Assert(((*it)->subdomain_id() == tria->locally_owned_subdomain() ||
86 tria->locally_owned_subdomain() == numbers::invalid_subdomain_id),
87 ExcInternalError());
88 }
89
91 cells[0]; // set iterator to master cell
92
93 // Store the slaves of each master
94 master2slaves[master_idx] = slaves;
95 // Save to which polygon this agglomerate correspond
96 master2polygon[master_idx] = n_agglomerations;
97
98 ++n_agglomerations; // an agglomeration has been performed, record it
99
100 create_bounding_box(cells); // fill the vector of bboxes
101
102 // Finally, return a polygonal iterator to the polytope just constructed.
103 return {cells[0], this};
104}
105
106template <int dim, int spacedim>
109 const AgglomerationContainer &cells,
110 const unsigned int fecollection_size)
111{
112 Assert(cells.size() > 0, ExcMessage("No cells to be agglomerated."));
113
114 if (cells.size() == 1)
115 hybrid_mesh = true; // mesh is made also by classical cells
116
117 // First index drives the selection of the master cell. After that, store the
118 // master cell.
119 const types::global_cell_index global_master_idx =
120 cells[0]->global_active_cell_index();
121 const types::global_cell_index master_idx = cells[0]->active_cell_index();
122 master_cells_container.push_back(cells[0]);
123 master_slave_relationships[global_master_idx] = -1;
124
125 const typename DoFHandler<dim>::active_cell_iterator cell_dh =
126 cells[0]->as_dof_handler_iterator(agglo_dh);
127 cell_dh->set_active_fe_index(CellAgglomerationType::master);
128
129 // Store slave cells and save the relationship with the parent
130 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>
131 slaves;
132 slaves.reserve(cells.size() - 1);
133 // exclude first cell since it's the master cell
134 for (auto it = ++cells.begin(); it != cells.end(); ++it)
135 {
136 slaves.push_back(*it);
137 master_slave_relationships[(*it)->global_active_cell_index()] =
138 global_master_idx; // mark each slave
139 master_slave_relationships_iterators[(*it)->active_cell_index()] =
140 cells[0];
141
142 const typename DoFHandler<dim>::active_cell_iterator cell =
143 (*it)->as_dof_handler_iterator(agglo_dh);
144 cell->set_active_fe_index(
145 fecollection_size); // slave cell (the last index)
146
147 // If we have a p::d::T, check that all cells are in the same subdomain.
148 // If serial, just check that the subdomain_id is invalid.
149 Assert(((*it)->subdomain_id() == tria->locally_owned_subdomain() ||
150 tria->locally_owned_subdomain() == numbers::invalid_subdomain_id),
151 ExcInternalError());
152 }
153
155 cells[0]; // set iterator to master cell
156
157 // Store the slaves of each master
158 master2slaves[master_idx] = slaves;
159 // Save to which polygon this agglomerate correspond
160 master2polygon[master_idx] = n_agglomerations;
161
162 ++n_agglomerations; // an agglomeration has been performed, record it
163
164 create_bounding_box(cells); // fill the vector of bboxes
165
166 // Finally, return a polygonal iterator to the polytope just constructed.
167 return {cells[0], this};
168}
169
170
171
172template <int dim, int spacedim>
173void
175 const AgglomerationContainer &cells)
176{
177 Assert(cells.size() > 0, ExcMessage("No cells to be agglomerated."));
178
179 // Compute the graph associated to the agglomerate
180 Utils::Graph g;
182
183 // The following vector will be filled with connected components
184 std::vector<std::vector<types::global_cell_index>> connected_components;
185 Utils::compute_connected_components(g, connected_components);
186
187 // Be verbose in debug mode
188#ifdef AGGLO_DEBUG
189 if (connected_components.size() > 1)
190 std::cout << "Disconnected elements. Connected components will be "
191 "computed and agglomerated together."
192 << std::endl;
193#endif
194
195 // Get connected components and define one agglomerate for each one of them.
196 std::vector<typename Triangulation<dim>::active_cell_iterator> agglomerate;
197 for (const auto &comp : connected_components)
198 {
199 agglomerate.reserve(comp.size());
200 for (const types::global_cell_index local_idx : comp)
201 agglomerate.push_back(cells[local_idx]);
202
203 // Perform agglomeration
204 define_agglomerate(agglomerate);
205 agglomerate.clear();
206 }
207}
208
209
210template <int dim, int spacedim>
211void
213 const Quadrature<dim> &cell_quadrature,
214 const UpdateFlags &flags,
215 const Quadrature<dim - 1> &face_quadrature,
216 const UpdateFlags &face_flags)
217{
218 agglomeration_quad = cell_quadrature;
219 agglomeration_flags = flags;
220 agglomeration_face_quad = face_quadrature;
222
223
224 no_values =
225 std::make_unique<FEValues<dim>>(*mapping,
226 dummy_fe,
229 update_JxW_values); // only for quadrature
230 no_face_values = std::make_unique<FEFaceValues<dim>>(
231 *mapping,
232 dummy_fe,
235 update_normal_vectors); // only for quadrature
236}
237
238template <int dim, int spacedim>
239void
241 const hp::QCollection<dim> &cell_qcollection,
242 const UpdateFlags &flags,
243 const hp::QCollection<dim - 1> &face_qcollection,
244 const UpdateFlags &face_flags)
245{
246 agglomeration_quad_collection = cell_qcollection;
247 agglomeration_flags = flags;
250
253 hp_no_values = std::make_unique<hp::FEValues<dim>>(
257 update_quadrature_points | update_JxW_values); // only for quadrature
258
259 hp_no_face_values = std::make_unique<hp::FEFaceValues<dim>>(
264 update_normal_vectors); // only for quadrature
265}
266
267
268
269template <int dim, int spacedim>
270unsigned int
274 unsigned int n_neighbors = 0;
275 for (const auto &f : cell->face_indices())
276 {
277 const auto &neighboring_cell = cell->neighbor(f);
278 if ((cell->face(f)->at_boundary()) ||
279 (neighboring_cell->is_active() &&
280 !are_cells_agglomerated(cell, neighboring_cell)))
281 {
282 ++n_neighbors;
283 }
284 }
285 return n_neighbors;
286}
287
288
290template <int dim, int spacedim>
291void
293 const std::unique_ptr<GridTools::Cache<dim, spacedim>> &cache_tria)
294{
295 tria = &cache_tria->get_triangulation();
296 mapping = &cache_tria->get_mapping();
297
298 agglo_dh.reinit(*tria);
299
300 if (const auto parallel_tria = dynamic_cast<
301 const dealii::parallel::TriangulationBase<dim, spacedim> *>(&*tria))
302 {
303 const std::weak_ptr<const Utilities::MPI::Partitioner> cells_partitioner =
304 parallel_tria->global_active_cell_index_partitioner();
306 cells_partitioner.lock()->locally_owned_range(), communicator);
307 }
308 else
309 {
310 master_slave_relationships.reinit(tria->n_active_cells(), MPI_COMM_SELF);
311 }
312
313 polytope_cache.clear();
314 bboxes.clear();
315
316 // First, update the pointer
317 cached_tria = std::make_unique<GridTools::Cache<dim, spacedim>>(
318 cache_tria->get_triangulation(), cache_tria->get_mapping());
319
322}
323
324
325
326template <int dim, int spacedim>
327void
329 const FiniteElement<dim> &fe_space)
330{
331 if (dynamic_cast<const FE_DGQ<dim> *>(&fe_space))
332 fe = std::make_unique<FE_DGQ<dim>>(fe_space.degree);
333 else if (dynamic_cast<const FE_AggloDGP<dim> *>(&fe_space))
334 fe = std::make_unique<FE_AggloDGP<dim>>(fe_space.degree);
335 else if (dynamic_cast<const FE_SimplexDGP<dim> *>(&fe_space))
336 fe = std::make_unique<FE_SimplexDGP<dim>>(fe_space.degree);
337 else
339 false,
340 ExcNotImplemented(
341 "Currently, this interface supports only DGQ and DGP bases."));
342
343 box_mapping = std::make_unique<MappingBox<dim>>(
344 bboxes,
345 master2polygon); // construct bounding box mapping
346
347 if (hybrid_mesh)
348 {
349 // the mesh is composed by standard and agglomerate cells. initialize
350 // classes needed for standard cells in order to treat that finite
351 // element space as defined on a standard shape and not on the
352 // BoundingBox.
354 std::make_unique<ScratchData>(*mapping,
355 *fe,
356 QGauss<dim>(2 * fe_space.degree + 2),
358 }
359
360
361 fe_collection.push_back(*fe); // master
362 fe_collection.push_back(
363 FE_Nothing<dim, spacedim>(fe->reference_cell())); // slave
364
366
367 // in case the tria is distributed, communicate ghost information with
368 // neighboring ranks
369 const bool needs_ghost_info =
370 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(&*tria) !=
371 nullptr;
372 if (needs_ghost_info)
374
376
377 if (needs_ghost_info)
379}
380
381template <int dim, int spacedim>
382void
384 const hp::FECollection<dim, spacedim> &fe_collection_in)
385{
386 is_hp_collection = true;
387
388 hp_fe_collection = std::make_unique<hp::FECollection<dim, spacedim>>(
389 fe_collection_in); // copy the input collection
390
391 box_mapping = std::make_unique<MappingBox<dim>>(
392 bboxes,
393 master2polygon); // construct bounding box mapping
394
395
396 if (hybrid_mesh)
397 {
398 AssertThrow(false,
399 ExcNotImplemented(
400 "Hybrid mesh is not implemented for hp::FECollection."));
401 }
402
403 for (unsigned int i = 0; i < fe_collection_in.size(); ++i)
404 {
405 if (dynamic_cast<const FESystem<dim> *>(&fe_collection_in[i]))
406 {
407 // System case
408 for (unsigned int b = 0; b < fe_collection_in[i].n_base_elements();
409 ++b)
410 {
411 if (!(dynamic_cast<const FE_DGQ<dim> *>(
412 &fe_collection_in[i].base_element(b)) ||
413 dynamic_cast<const FE_AggloDGP<dim> *>(
414 &fe_collection_in[i].base_element(b)) ||
415 dynamic_cast<const FE_SimplexDGP<dim> *>(
416 &fe_collection_in[i].base_element(b)) ||
417 dynamic_cast<const FE_Nothing<dim> *>(
418 &fe_collection_in[i].base_element(b))))
420 false,
421 ExcNotImplemented(
422 "Currently, this interface supports only DGQ and DGP bases."));
423 }
425 else
426 {
427 // Scalar case
428 if (!(dynamic_cast<const FE_DGQ<dim> *>(&fe_collection_in[i]) ||
429 dynamic_cast<const FE_AggloDGP<dim> *>(&fe_collection_in[i]) ||
430 dynamic_cast<const FE_SimplexDGP<dim> *>(&fe_collection_in[i])))
432 false,
433 ExcNotImplemented(
434 "Currently, this interface supports only DGQ and DGP bases."));
435 }
436 fe_collection.push_back(fe_collection_in[i]);
437 }
438
439 Assert(fe_collection[0].n_components() >= 1,
440 ExcMessage("Invalid FE: must have at least one component."));
441 if (fe_collection[0].n_components() == 1)
442 {
444 }
445 else if (fe_collection[0].n_components() > 1)
446 {
447 std::vector<const FiniteElement<dim, spacedim> *> base_elements;
448 std::vector<unsigned int> multiplicities;
449 for (unsigned int b = 0; b < fe_collection[0].n_base_elements(); ++b)
450 {
451 base_elements.push_back(new FE_Nothing<dim, spacedim>());
452 multiplicities.push_back(fe_collection[0].element_multiplicity(b));
453 }
454 FESystem<dim, spacedim> fe_system_nothing(base_elements, multiplicities);
455 for (const auto *ptr : base_elements)
456 delete ptr;
457 fe_collection.push_back(fe_system_nothing);
458 }
459
461
462 // in case the tria is distributed, communicate ghost information with
463 // neighboring ranks
464 const bool needs_ghost_info =
465 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(&*tria) !=
466 nullptr;
467 if (needs_ghost_info)
469
471
472 if (needs_ghost_info)
474}
475
476template <int dim, int spacedim>
477void
479 const AgglomerationContainer &polytope)
480{
482 ExcMessage("No agglomeration has been performed."));
483 Assert(dim > 1, ExcNotImplemented());
484
485 std::vector<Point<spacedim>> pts; // store all the vertices
486 for (const auto &cell : polytope)
487 for (const auto i : cell->vertex_indices())
488 pts.push_back(cell->vertex(i));
489
490 bboxes.emplace_back(pts);
491}
492
493
494
495template <int dim, int spacedim>
496void
498{
500 ExcMessage("No agglomeration has been performed."));
501 Assert(
502 agglo_dh.n_dofs() > 0,
503 ExcMessage(
504 "The DoFHandler associated to the agglomeration has not been initialized."
505 "It's likely that you forgot to distribute the DoFs. You may want"
506 "to check if a call to `initialize_hp_structure()` has been done."));
507
509 for (const auto &cell : master_cells_container)
510 {
513 }
514
516 {
517 // communicate the number of faces
520 // send information about boundaries and neighboring polytopes id
523
526 }
527}
528
529
530
531template <int dim, int spacedim>
532void
534{
535 const unsigned int dofs_per_cell = fe->dofs_per_cell;
536 for (const auto &polytope : polytope_iterators())
537 {
538 if (polytope->is_locally_owned())
539 {
540 const unsigned int n_faces = polytope->n_faces();
541 for (unsigned int f = 0; f < n_faces; ++f)
542 {
543 if (!polytope->at_boundary(f))
544 {
545 const auto &neigh_polytope = polytope->neighbor(f);
546 if (!neigh_polytope->is_locally_owned())
547 {
548 // Neighboring polytope is ghosted.
549
550 // Compute shape functions at the interface
551 const auto &current_fe = reinit(polytope, f);
552
553 std::vector<Point<spacedim>> qpoints_to_send =
554 current_fe.get_quadrature_points();
555
556 const std::vector<double> &jxws_to_send =
557 current_fe.get_JxW_values();
558
559 const std::vector<Tensor<1, spacedim>> &normals_to_send =
560 current_fe.get_normal_vectors();
561
562
563 const types::subdomain_id neigh_rank =
564 neigh_polytope->subdomain_id();
565
566 std::pair<CellId, unsigned int> cell_and_face{
567 polytope->id(), f};
568 // Prepare data to send
569 local_qpoints[neigh_rank].emplace(cell_and_face,
570 qpoints_to_send);
571
572 local_jxws[neigh_rank].emplace(cell_and_face,
573 jxws_to_send);
574
575 local_normals[neigh_rank].emplace(cell_and_face,
576 normals_to_send);
577
578
579 const unsigned int n_qpoints = qpoints_to_send.size();
580
581 // TODO: check `agglomeration_flags` before computing
582 // values and gradients.
583 std::vector<std::vector<double>> values_per_qpoints(
584 dofs_per_cell);
585
586 std::vector<std::vector<Tensor<1, spacedim>>>
587 gradients_per_qpoints(dofs_per_cell);
588
589 for (unsigned int i = 0; i < dofs_per_cell; ++i)
590 {
591 values_per_qpoints[i].resize(n_qpoints);
592 gradients_per_qpoints[i].resize(n_qpoints);
593 for (unsigned int q = 0; q < n_qpoints; ++q)
594 {
595 values_per_qpoints[i][q] =
596 current_fe.shape_value(i, q);
597 gradients_per_qpoints[i][q] =
598 current_fe.shape_grad(i, q);
599 }
600 }
601
602 local_values[neigh_rank].emplace(cell_and_face,
603 values_per_qpoints);
604 local_gradients[neigh_rank].emplace(
605 cell_and_face, gradients_per_qpoints);
606 }
607 }
608 }
609 }
610 }
611
612 // Finally, exchange with neighboring ranks
613 recv_qpoints = Utilities::MPI::some_to_some(communicator, local_qpoints);
614 recv_jxws = Utilities::MPI::some_to_some(communicator, local_jxws);
615 recv_normals = Utilities::MPI::some_to_some(communicator, local_normals);
616 recv_values = Utilities::MPI::some_to_some(communicator, local_values);
617 recv_gradients = Utilities::MPI::some_to_some(communicator, local_gradients);
619
620
621
622template <int dim, int spacedim>
626 &cells,
628 &master_cell) const
629{
630 Assert(is_master_cell(master_cell),
631 ExcMessage("This must be a master cell."));
632
633 std::vector<Point<dim>> vec_pts;
634 std::vector<double> vec_JxWs;
635
636 if (!is_hp_collection)
637 {
638 // Original version: handle case without hp::FECollection
639 for (const auto &dummy_cell : cells)
640 {
641 no_values->reinit(dummy_cell);
642 auto q_points = no_values->get_quadrature_points(); // real qpoints
643 const auto &JxWs = no_values->get_JxW_values();
645 std::transform(q_points.begin(),
646 q_points.end(),
647 std::back_inserter(vec_pts),
648 [&](const Point<spacedim> &p) { return p; });
649 std::transform(JxWs.begin(),
650 JxWs.end(),
651 std::back_inserter(vec_JxWs),
652 [&](const double w) { return w; });
653 }
654 }
655 else
656 {
657 // Handle the hp::FECollection case
658 const auto &master_cell_as_dh_iterator =
659 master_cell->as_dof_handler_iterator(agglo_dh);
660 for (const auto &dummy_cell : cells)
661 {
662 // The following verbose call is necessary to handle cases where
663 // different slave cells on different polytopes use different
664 // quadrature rules. If the hp::QCollection contains multiple
665 // elements, calling hp_no_values->reinit(dummy_cell) won't work
666 // because it cannot infer the correct quadrature rule. By explicitly
667 // passing the active FE index as q_index, and setting mapping_index
668 // and fe_index to 0, we ensure that the dummy cell uses the same
669 // quadrature rule as its corresponding master cell. This assumes a
670 // one-to-one correspondence between hp::QCollection and
671 // hp::FECollection, which is the convention in deal.II. However, this
672 // implementation does not support cases where hp::QCollection and
673 // hp::FECollection have different sizes.
674 // TODO: Refactor the architecture to better handle numerical
675 // integration for hp::QCollection.
676 hp_no_values->reinit(dummy_cell,
677 master_cell_as_dh_iterator->active_fe_index(),
678 0,
679 0);
680 auto q_points = hp_no_values->get_present_fe_values()
681 .get_quadrature_points(); // real qpoints
682 const auto &JxWs =
683 hp_no_values->get_present_fe_values().get_JxW_values();
684
685 std::transform(q_points.begin(),
686 q_points.end(),
687 std::back_inserter(vec_pts),
688 [&](const Point<spacedim> &p) { return p; });
689 std::transform(JxWs.begin(),
690 JxWs.end(),
691 std::back_inserter(vec_JxWs),
692 [&](const double w) { return w; });
693 }
694 }
695
696 // Map back each point in real space by using the map associated to the
697 // bounding box.
698 std::vector<Point<dim>> unit_points(vec_pts.size());
699 const auto &bbox =
700 bboxes[master2polygon.at(master_cell->active_cell_index())];
701 unit_points.reserve(vec_pts.size());
702
703 for (unsigned int i = 0; i < vec_pts.size(); i++)
704 unit_points[i] = bbox.real_to_unit(vec_pts[i]);
705
706 return Quadrature<dim>(unit_points, vec_JxWs);
707}
708
709
710
711template <int dim, int spacedim>
712void
714{
715 Assert(agglo_dh.get_triangulation().n_cells() > 0,
716 ExcMessage(
717 "Triangulation must not be empty upon calling this function."));
719 ExcMessage("No agglomeration has been performed."));
720
721 agglo_dh.distribute_dofs(fe_collection);
722 // euler_mapping = std::make_unique<
723 // MappingFEField<dim, spacedim,
724 // LinearAlgebra::distributed::Vector<double>>>( euler_dh, euler_vector);
725}
726
727
728
729template <int dim, int spacedim>
732 const AgglomerationIterator<dim, spacedim> &polytope) const
733{
734 // Assert(euler_mapping,
735 // ExcMessage("The mapping describing the physical element stemming
736 // from "
737 // "agglomeration has not been set up."));
738
739 const auto &deal_cell = polytope->as_dof_handler_iterator(agglo_dh);
740
741 // First check if the polytope is made just by a single cell. If so, use
742 // classical FEValues
743 // if (polytope->n_background_cells() == 1)
744 // return standard_scratch->reinit(deal_cell);
745
746 const auto &agglo_cells = polytope->get_agglomerate();
747
748 Quadrature<dim> agglo_quad = agglomerated_quadrature(agglo_cells, deal_cell);
749
750 if (!is_hp_collection)
751 {
752 // Original version: handle case without hp::FECollection
753 agglomerated_scratch = std::make_unique<ScratchData>(*box_mapping,
754 fe_collection[0],
755 agglo_quad,
757 }
758 else
759 {
760 // Handle the hp::FECollection case
761 agglomerated_scratch = std::make_unique<ScratchData>(*box_mapping,
762 polytope->get_fe(),
763 agglo_quad,
765 }
766 return agglomerated_scratch->reinit(deal_cell);
767}
768
769
770
771template <int dim, int spacedim>
775 const unsigned int face_index,
777 &agglo_isv_ptr) const
778{
780 reinit_master(cell, face_index, agglo_isv_ptr, *this);
781}
782
783
784
785template <int dim, int spacedim>
789 const unsigned int face_index) const
790{
791 // Assert(euler_mapping,
792 // ExcMessage("The mapping describing the physical element stemming
793 // from "
794 // "agglomeration has not been set up."));
795
796 const auto &deal_cell = polytope->as_dof_handler_iterator(agglo_dh);
797 Assert(is_master_cell(deal_cell), ExcMessage("This should be true."));
798
800 reinit_master(deal_cell, face_index, agglomerated_isv_bdary, *this);
801}
802
803
804
805template <int dim, int spacedim>
806std::pair<const FEValuesBase<dim, spacedim> &,
809 const AgglomerationIterator<dim, spacedim> &polytope_in,
810 const AgglomerationIterator<dim, spacedim> &neigh_polytope,
811 const unsigned int local_in,
812 const unsigned int local_neigh) const
813{
814 // If current and neighboring polytopes are both locally owned, then compute
815 // the jump in the classical way without needing information about ghosted
816 // entities.
817 if (polytope_in->is_locally_owned() && neigh_polytope->is_locally_owned())
818 {
819 const auto &cell_in = polytope_in->as_dof_handler_iterator(agglo_dh);
820 const auto &neigh_cell =
821 neigh_polytope->as_dof_handler_iterator(agglo_dh);
822
823 const auto &fe_in =
825 reinit_master(cell_in, local_in, agglomerated_isv, *this);
826 const auto &fe_out =
828 reinit_master(neigh_cell, local_neigh, agglomerated_isv_neigh, *this);
829 std::pair<const FEValuesBase<dim, spacedim> &,
831 my_p(fe_in, fe_out);
832
833 return my_p;
834 }
835 else
836 {
837 Assert((polytope_in->is_locally_owned() &&
838 !neigh_polytope->is_locally_owned()),
839 ExcInternalError());
840
841 const auto &cell = polytope_in->as_dof_handler_iterator(agglo_dh);
842 const auto &bbox = bboxes[master2polygon.at(cell->active_cell_index())];
843 // const double bbox_measure = bbox.volume();
844
845 const unsigned int neigh_rank = neigh_polytope->subdomain_id();
846 const CellId &neigh_id = neigh_polytope->id();
847
848 // Retrieve qpoints,JxWs, normals sent previously from the neighboring
849 // rank.
850 std::vector<Point<spacedim>> &real_qpoints =
851 recv_qpoints.at(neigh_rank).at({neigh_id, local_neigh});
852
853 const auto &JxWs = recv_jxws.at(neigh_rank).at({neigh_id, local_neigh});
854
855 std::vector<Tensor<1, spacedim>> &normals =
856 recv_normals.at(neigh_rank).at({neigh_id, local_neigh});
857
858 // Apply the necessary scalings due to the bbox.
859 std::vector<Point<spacedim>> final_unit_q_points;
860 std::transform(real_qpoints.begin(),
861 real_qpoints.end(),
862 std::back_inserter(final_unit_q_points),
863 [&](const Point<spacedim> &p) {
864 return bbox.real_to_unit(p);
865 });
866
867 // std::vector<double> scale_factors(final_unit_q_points.size());
868 // std::vector<double> scaled_weights(final_unit_q_points.size());
869 // std::vector<Tensor<1, dim>> scaled_normals(final_unit_q_points.size());
870
871 // Since we received normal vectors from a neighbor, we have to swap
872 // the
873 // // sign of the vector in order to have outward normals.
874 // for (unsigned int q = 0; q < final_unit_q_points.size(); ++q)
875 // {
876 // for (unsigned int direction = 0; direction < spacedim; ++direction)
877 // scaled_normals[q][direction] =
878 // normals[q][direction] * (bbox.side_length(direction));
879
880 // scaled_normals[q] *= -1;
881
882 // scaled_weights[q] =
883 // (JxWs[q] * scaled_normals[q].norm()) / bbox_measure;
884 // scaled_normals[q] /= scaled_normals[q].norm();
885 // }
886 for (unsigned int q = 0; q < final_unit_q_points.size(); ++q)
887 normals[q] *= -1;
888
889
891 final_unit_q_points, JxWs, normals);
892
894 std::make_unique<NonMatching::FEImmersedSurfaceValues<spacedim>>(
895 *box_mapping, *fe, surface_quad, agglomeration_face_flags);
896
897
898 agglomerated_isv->reinit(cell);
899
900 std::pair<const FEValuesBase<dim, spacedim> &,
903
904 return my_p;
905 }
906}
907
908
909
910template <int dim, int spacedim>
911template <typename SparsityPatternType, typename Number>
912void
914 SparsityPatternType &dsp,
915 const AffineConstraints<Number> &constraints,
916 const bool keep_constrained_dofs,
917 const types::subdomain_id subdomain_id)
918{
920 ExcMessage("The agglomeration has not been set up correctly."));
921 Assert(dsp.empty(),
922 ExcMessage(
923 "The Sparsity pattern must be empty upon calling this function."));
924
925 const IndexSet &locally_owned_dofs = agglo_dh.locally_owned_dofs();
926 const IndexSet locally_relevant_dofs =
928
929 if constexpr (std::is_same_v<SparsityPatternType, DynamicSparsityPattern>)
930 dsp.reinit(locally_owned_dofs.size(),
931 locally_owned_dofs.size(),
932 locally_relevant_dofs);
933 else if constexpr (std::is_same_v<SparsityPatternType,
935 dsp.reinit(locally_owned_dofs, communicator);
936 else
937 AssertThrow(false, ExcNotImplemented());
938
939 // Create the sparsity pattern corresponding only to volumetric terms. The
940 // fluxes needed by DG methods will be filled later.
942 agglo_dh, dsp, constraints, keep_constrained_dofs, subdomain_id);
943
944
945 if (!is_hp_collection)
946 {
947 // Original version: handle case without hp::FECollection
948 const unsigned int dofs_per_cell = agglo_dh.get_fe(0).n_dofs_per_cell();
949 std::vector<types::global_dof_index> current_dof_indices(dofs_per_cell);
950 std::vector<types::global_dof_index> neighbor_dof_indices(dofs_per_cell);
951
952 // Loop over all locally owned polytopes, find the neighbor (also ghosted)
953 // and add fluxes to the sparsity pattern.
954 for (const auto &polytope : polytope_iterators())
955 {
956 if (polytope->is_locally_owned())
957 {
958 const unsigned int n_current_faces = polytope->n_faces();
959 polytope->get_dof_indices(current_dof_indices);
960 for (unsigned int f = 0; f < n_current_faces; ++f)
961 {
962 const auto &neigh_polytope = polytope->neighbor(f);
963 if (neigh_polytope.state() == IteratorState::valid)
964 {
965 neigh_polytope->get_dof_indices(neighbor_dof_indices);
966 constraints.add_entries_local_to_global(
967 current_dof_indices,
968 neighbor_dof_indices,
969 dsp,
970 keep_constrained_dofs,
971 {});
972 }
973 }
974 }
975 }
976 }
977 else
978 {
979 // Handle the hp::FECollection case
980
981 // Loop over all locally owned polytopes, find the neighbor (also ghosted)
982 // and add fluxes to the sparsity pattern.
983 for (const auto &polytope : polytope_iterators())
984 {
985 if (polytope->is_locally_owned())
986 {
987 const unsigned int current_dofs_per_cell =
988 polytope->get_fe().dofs_per_cell;
989 std::vector<types::global_dof_index> current_dof_indices(
990 current_dofs_per_cell);
991
992 const unsigned int n_current_faces = polytope->n_faces();
993 polytope->get_dof_indices(current_dof_indices);
994 for (unsigned int f = 0; f < n_current_faces; ++f)
995 {
996 const auto &neigh_polytope = polytope->neighbor(f);
997 if (neigh_polytope.state() == IteratorState::valid)
998 {
999 const unsigned int neighbor_dofs_per_cell =
1000 neigh_polytope->get_fe().dofs_per_cell;
1001 std::vector<types::global_dof_index> neighbor_dof_indices(
1002 neighbor_dofs_per_cell);
1003
1004 neigh_polytope->get_dof_indices(neighbor_dof_indices);
1005 constraints.add_entries_local_to_global(
1006 current_dof_indices,
1007 neighbor_dof_indices,
1008 dsp,
1009 keep_constrained_dofs,
1010 {});
1011 }
1012 }
1013 }
1014 }
1015 }
1016
1017
1018
1019 if constexpr (std::is_same_v<SparsityPatternType,
1021 dsp.compress();
1022}
1023
1024
1025
1026template <int dim, int spacedim>
1027void
1029{
1030 [[maybe_unused]] const auto parallel_triangulation =
1031 dynamic_cast<const parallel::TriangulationBase<dim, spacedim> *>(&*tria);
1032 Assert(parallel_triangulation != nullptr, ExcInternalError());
1033
1034 const unsigned int n_dofs_per_cell = fe->dofs_per_cell;
1035 std::vector<types::global_dof_index> global_dof_indices(n_dofs_per_cell);
1036 for (const auto &polytope : polytope_iterators())
1037 if (polytope->is_locally_owned())
1038 {
1039 const CellId &master_cell_id = polytope->id();
1040
1041 const auto polytope_dh = polytope->as_dof_handler_iterator(agglo_dh);
1042 polytope_dh->get_dof_indices(global_dof_indices);
1043
1044
1045 const auto &agglomerate = polytope->get_agglomerate();
1046
1047 for (const auto &cell : agglomerate)
1048 {
1049 // interior, locally owned, cell
1050 for (const auto &f : cell->face_indices())
1051 {
1052 if (!cell->at_boundary(f))
1053 {
1054 const auto &neighbor = cell->neighbor(f);
1055 if (neighbor->is_ghost())
1056 {
1057 // key of the map: the rank to which send the data
1058 const types::subdomain_id neigh_rank =
1059 neighbor->subdomain_id();
1060
1061 // inform the "standard" neighbor about the neighboring
1062 // id and its master cell
1063 local_cell_ids_neigh_cell[neigh_rank].emplace(
1064 cell->id(), master_cell_id);
1065
1066 // inform the neighboring rank that this master cell
1067 // (hence polytope) has the following DoF indices
1068 local_ghost_dofs[neigh_rank].emplace(
1069 master_cell_id, global_dof_indices);
1070
1071 // ...same for bounding boxes
1072 const auto &bbox = bboxes[polytope->index()];
1073 local_ghosted_bbox[neigh_rank].emplace(master_cell_id,
1074 bbox);
1075 }
1076 }
1077 }
1078 }
1079 }
1080
1083
1084 // Exchange with neighboring ranks the neighboring bounding boxes
1087
1088 // Exchange with neighboring ranks the neighboring ghosted DoFs
1091}
1092
1093
1094
1095namespace dealii
1096{
1097 namespace internal
1098 {
1099 template <int dim, int spacedim>
1101 {
1102 public:
1103 static const FEValuesBase<dim, spacedim> &
1106 const unsigned int face_index,
1108 &agglo_isv_ptr,
1110 {
1111 Assert(handler.is_master_cell(cell),
1112 ExcMessage("This cell must be a master one."));
1113
1114 AgglomerationIterator<dim, spacedim> it{cell, &handler};
1115 const auto &neigh_polytope = it->neighbor(face_index);
1116
1117 const CellId polytope_in_id = cell->id();
1118
1119 // Retrieve the bounding box of the agglomeration
1120 const auto &bbox =
1121 handler.bboxes[handler.master2polygon.at(cell->active_cell_index())];
1122
1123 CellId polytope_out_id;
1124 if (neigh_polytope.state() == IteratorState::valid)
1125 polytope_out_id = neigh_polytope->id();
1126 else
1127 polytope_out_id = polytope_in_id; // on the boundary. Same id
1128
1129 const auto &common_face = handler.polytope_cache.interface.at(
1130 {polytope_in_id, polytope_out_id});
1131
1132 std::vector<Point<spacedim>> final_unit_q_points;
1133 std::vector<double> final_weights;
1134 std::vector<Tensor<1, dim>> final_normals;
1135
1136 if (!handler.is_hp_collection)
1137 {
1138 // Original version: handle case without hp::FECollection
1139 const unsigned int expected_qpoints =
1140 common_face.size() * handler.agglomeration_face_quad.size();
1141 final_unit_q_points.reserve(expected_qpoints);
1142 final_weights.reserve(expected_qpoints);
1143 final_normals.reserve(expected_qpoints);
1144
1145
1146 for (const auto &[deal_cell, local_face_idx] : common_face)
1147 {
1148 handler.no_face_values->reinit(deal_cell, local_face_idx);
1149
1150 const auto &q_points =
1151 handler.no_face_values->get_quadrature_points();
1152 const auto &JxWs = handler.no_face_values->get_JxW_values();
1153 const auto &normals =
1154 handler.no_face_values->get_normal_vectors();
1155
1156 const unsigned int n_qpoints_agglo = q_points.size();
1157
1158 for (unsigned int q = 0; q < n_qpoints_agglo; ++q)
1159 {
1160 final_unit_q_points.push_back(
1161 bbox.real_to_unit(q_points[q]));
1162 final_weights.push_back(JxWs[q]);
1163 final_normals.push_back(normals[q]);
1164 }
1165 }
1166 }
1167 else
1168 {
1169 // Handle the hp::FECollection case
1170 unsigned int higher_order_quad_index = cell->active_fe_index();
1171 if (neigh_polytope.state() == IteratorState::valid)
1172 if (handler
1173 .agglomeration_face_quad_collection[cell->active_fe_index()]
1174 .size() <
1175 handler
1177 ->active_fe_index()]
1178 .size())
1179 higher_order_quad_index = neigh_polytope->active_fe_index();
1180
1181 const unsigned int expected_qpoints =
1182 common_face.size() *
1183 handler
1184 .agglomeration_face_quad_collection[higher_order_quad_index]
1185 .size();
1186 final_unit_q_points.reserve(expected_qpoints);
1187 final_weights.reserve(expected_qpoints);
1188 final_normals.reserve(expected_qpoints);
1189
1190 for (const auto &[deal_cell, local_face_idx] : common_face)
1191 {
1192 handler.hp_no_face_values->reinit(
1193 deal_cell, local_face_idx, higher_order_quad_index, 0, 0);
1194
1195 const auto &q_points =
1196 handler.hp_no_face_values->get_present_fe_values()
1197 .get_quadrature_points();
1198 const auto &JxWs =
1199 handler.hp_no_face_values->get_present_fe_values()
1200 .get_JxW_values();
1201 const auto &normals =
1202 handler.hp_no_face_values->get_present_fe_values()
1203 .get_normal_vectors();
1204
1205 const unsigned int n_qpoints_agglo = q_points.size();
1206
1207 for (unsigned int q = 0; q < n_qpoints_agglo; ++q)
1208 {
1209 final_unit_q_points.push_back(
1210 bbox.real_to_unit(q_points[q]));
1211 final_weights.push_back(JxWs[q]);
1212 final_normals.push_back(normals[q]);
1213 }
1214 }
1215 }
1216
1217
1219 final_unit_q_points, final_weights, final_normals);
1220
1221 if (!handler.is_hp_collection)
1222 {
1223 agglo_isv_ptr =
1224 std::make_unique<NonMatching::FEImmersedSurfaceValues<spacedim>>(
1225 *(handler.box_mapping),
1226 *(handler.fe),
1227 surface_quad,
1228 handler.agglomeration_face_flags);
1229 }
1230 else
1231 {
1232 agglo_isv_ptr =
1233 std::make_unique<NonMatching::FEImmersedSurfaceValues<spacedim>>(
1234 *(handler.box_mapping),
1235 cell->get_fe(),
1236 surface_quad,
1237 handler.agglomeration_face_flags);
1238 }
1239
1240 agglo_isv_ptr->reinit(cell);
1241
1242 return *agglo_isv_ptr;
1243 }
1244
1245
1246
1253 static void
1256 &master_cell,
1258 {
1259 Assert(
1260 handler.master_slave_relationships[master_cell
1261 ->global_active_cell_index()] ==
1262 -1,
1263 ExcMessage("The present cell with index " +
1264 std::to_string(master_cell->global_active_cell_index()) +
1265 "is not a master one."));
1266
1267 const auto &agglomeration = handler.get_agglomerate(master_cell);
1268 const types::global_cell_index current_polytope_index =
1269 handler.master2polygon.at(master_cell->active_cell_index());
1270
1271 CellId current_polytope_id = master_cell->id();
1272
1273
1274 std::set<types::global_cell_index> visited_polygonal_neighbors;
1275
1276 std::map<unsigned int, CellId> face_to_neigh_id;
1277
1278 std::map<unsigned int, bool> is_face_at_boundary;
1279
1280 // same as above, but with CellId
1281 std::set<CellId> visited_polygonal_neighbors_id;
1282 unsigned int ghost_counter = 0;
1283
1284 for (const auto &cell : agglomeration)
1285 {
1287 cell->active_cell_index();
1288
1289 const CellId cell_id = cell->id();
1290
1291 for (const auto f : cell->face_indices())
1292 {
1293 const auto &neighboring_cell = cell->neighbor(f);
1294
1295 const bool valid_neighbor =
1296 neighboring_cell.state() == IteratorState::valid;
1297
1298 if (valid_neighbor)
1299 {
1300 if (neighboring_cell->is_locally_owned() &&
1301 !handler.are_cells_agglomerated(cell, neighboring_cell))
1302 {
1303 // - cell is not on the boundary,
1304 // - it's not agglomerated with the neighbor. If so,
1305 // it's a neighbor of the present agglomeration
1306 // std::cout << " (from rank) "
1307 // << Utilities::MPI::this_mpi_process(
1308 // handler.communicator)
1309 // << std::endl;
1310
1311 // std::cout
1312 // << "neighbor locally owned? " << std::boolalpha
1313 // << neighboring_cell->is_locally_owned() <<
1314 // std::endl;
1315 // if (neighboring_cell->is_ghost())
1316 // handler.ghosted_indices.push_back(
1317 // neighboring_cell->active_cell_index());
1318
1319 // a new face of the agglomeration has been
1320 // discovered.
1321 handler.polygon_boundary[master_cell].push_back(
1322 cell->face(f));
1323
1324 // global index of neighboring deal.II cell
1325 const types::global_cell_index neighboring_cell_index =
1326 neighboring_cell->active_cell_index();
1327
1328 // master cell for the neighboring polytope
1329 const auto &master_of_neighbor =
1331 neighboring_cell_index);
1332
1333 const auto nof = cell->neighbor_of_neighbor(f);
1334
1335 if (handler.is_slave_cell(neighboring_cell))
1336 {
1337 // index of the neighboring polytope
1339 neighbor_polytope_index =
1340 handler.master2polygon.at(
1341 master_of_neighbor->active_cell_index());
1342
1343 CellId neighbor_polytope_id =
1344 master_of_neighbor->id();
1345
1346 if (visited_polygonal_neighbors.find(
1347 neighbor_polytope_index) ==
1348 std::end(visited_polygonal_neighbors))
1349 {
1350 // found a neighbor
1351
1352 const unsigned int n_face =
1354 [current_polytope_index];
1355
1356 handler.polytope_cache.cell_face_at_boundary[{
1357 current_polytope_index, n_face}] = {
1358 false, master_of_neighbor};
1359
1360 is_face_at_boundary[n_face] = true;
1361
1363 [current_polytope_index];
1364
1365 visited_polygonal_neighbors.insert(
1366 neighbor_polytope_index);
1367 }
1368
1369
1370 if (handler.polytope_cache.visited_cell_and_faces
1371 .find({cell_index, f}) ==
1372 std::end(handler.polytope_cache
1373 .visited_cell_and_faces))
1374 {
1375 handler.polytope_cache
1376 .interface[{
1377 current_polytope_id, neighbor_polytope_id}]
1378 .emplace_back(cell, f);
1379
1380 handler.polytope_cache.visited_cell_and_faces
1381 .insert({cell_index, f});
1382 }
1383
1384
1385 if (handler.polytope_cache.visited_cell_and_faces
1386 .find({neighboring_cell_index, nof}) ==
1387 std::end(handler.polytope_cache
1388 .visited_cell_and_faces))
1389 {
1390 handler.polytope_cache
1391 .interface[{
1392 neighbor_polytope_id, current_polytope_id}]
1393 .emplace_back(neighboring_cell, nof);
1394
1395 handler.polytope_cache.visited_cell_and_faces
1396 .insert({neighboring_cell_index, nof});
1397 }
1398 }
1399 else
1400 {
1401 // neighboring cell is a master
1402
1403 // save the pair of neighboring cells
1405 neighbor_polytope_index =
1406 handler.master2polygon.at(
1407 neighboring_cell_index);
1408
1409 CellId neighbor_polytope_id =
1410 neighboring_cell->id();
1411
1412 if (visited_polygonal_neighbors.find(
1413 neighbor_polytope_index) ==
1414 std::end(visited_polygonal_neighbors))
1415 {
1416 // found a neighbor
1417 const unsigned int n_face =
1419 [current_polytope_index];
1420
1421
1422 handler.polytope_cache.cell_face_at_boundary[{
1423 current_polytope_index, n_face}] = {
1424 false, neighboring_cell};
1425
1426 is_face_at_boundary[n_face] = true;
1427
1429 [current_polytope_index];
1430
1431 visited_polygonal_neighbors.insert(
1432 neighbor_polytope_index);
1433 }
1434
1435
1436
1437 if (handler.polytope_cache.visited_cell_and_faces
1438 .find({cell_index, f}) ==
1439 std::end(handler.polytope_cache
1440 .visited_cell_and_faces))
1441 {
1442 handler.polytope_cache
1443 .interface[{
1444 current_polytope_id, neighbor_polytope_id}]
1445 .emplace_back(cell, f);
1446
1447 handler.polytope_cache.visited_cell_and_faces
1448 .insert({cell_index, f});
1449 }
1450
1451 if (handler.polytope_cache.visited_cell_and_faces
1452 .find({neighboring_cell_index, nof}) ==
1453 std::end(handler.polytope_cache
1454 .visited_cell_and_faces))
1455 {
1456 handler.polytope_cache
1457 .interface[{
1458 neighbor_polytope_id, current_polytope_id}]
1459 .emplace_back(neighboring_cell, nof);
1460
1461 handler.polytope_cache.visited_cell_and_faces
1462 .insert({neighboring_cell_index, nof});
1463 }
1464 }
1465 }
1466 else if (neighboring_cell->is_ghost())
1467 {
1468 const auto nof = cell->neighbor_of_neighbor(f);
1469
1470 // from neighboring rank,receive the association
1471 // between standard cell ids and neighboring polytope.
1472 // This tells to the current rank that the
1473 // neighboring cell has the following CellId as master
1474 // cell.
1475 const auto &check_neigh_poly_ids =
1476 handler.recv_cell_ids_neigh_cell.at(
1477 neighboring_cell->subdomain_id());
1478
1479 const CellId neighboring_cell_id =
1480 neighboring_cell->id();
1481
1482 const CellId &check_neigh_polytope_id =
1483 check_neigh_poly_ids.at(neighboring_cell_id);
1484
1485 // const auto master_index =
1486 // master_indices[ghost_counter];
1487
1488 if (visited_polygonal_neighbors_id.find(
1489 check_neigh_polytope_id) ==
1490 std::end(visited_polygonal_neighbors_id))
1491 {
1492 handler.polytope_cache.cell_face_at_boundary[{
1493 current_polytope_index,
1495 [current_polytope_index]}] = {false,
1496 neighboring_cell};
1497
1498
1499 // record the cell id of the neighboring polytope
1500 handler.polytope_cache.ghosted_master_id[{
1501 current_polytope_id,
1503 [current_polytope_index]}] =
1504 check_neigh_polytope_id;
1505
1506
1507 const unsigned int n_face =
1509 [current_polytope_index];
1510
1511 face_to_neigh_id[n_face] = check_neigh_polytope_id;
1512
1513 is_face_at_boundary[n_face] = false;
1514
1515
1516 // increment number of faces
1518 [current_polytope_index];
1519
1520 visited_polygonal_neighbors_id.insert(
1521 check_neigh_polytope_id);
1522
1523 // ghosted polytope has been found, increment
1524 // ghost counter
1525 ++ghost_counter;
1526 }
1527
1528
1529
1530 if (handler.polytope_cache.visited_cell_and_faces_id
1531 .find({cell_id, f}) ==
1532 std::end(
1533 handler.polytope_cache.visited_cell_and_faces_id))
1534 {
1535 handler.polytope_cache
1536 .interface[{
1537 current_polytope_id, check_neigh_polytope_id}]
1538 .emplace_back(cell, f);
1539
1540 // std::cout << "ADDED ("
1541 // << cell->active_cell_index() << ")
1542 // BETWEEN "
1543 // << current_polytope_id << " e "
1544 // << check_neigh_polytope_id <<
1545 // std::endl;
1546
1547 handler.polytope_cache.visited_cell_and_faces_id
1548 .insert({cell_id, f});
1549 }
1550
1551
1552 if (handler.polytope_cache.visited_cell_and_faces_id
1553 .find({neighboring_cell_id, nof}) ==
1554 std::end(
1555 handler.polytope_cache.visited_cell_and_faces_id))
1556 {
1557 handler.polytope_cache
1558 .interface[{
1559 check_neigh_polytope_id, current_polytope_id}]
1560 .emplace_back(neighboring_cell, nof);
1561
1562 handler.polytope_cache.visited_cell_and_faces_id
1563 .insert({neighboring_cell_id, nof});
1564 }
1565 }
1566 }
1567 else if (cell->face(f)->at_boundary())
1568 {
1569 // Boundary face of a boundary cell.
1570 // Note that the neighboring cell must be invalid.
1571
1572 handler.polygon_boundary[master_cell].push_back(
1573 cell->face(f));
1574
1575 if (visited_polygonal_neighbors.find(
1576 std::numeric_limits<unsigned int>::max()) ==
1577 std::end(visited_polygonal_neighbors))
1578 {
1579 // boundary face. Notice that `neighboring_cell` is
1580 // invalid here.
1581 handler.polytope_cache.cell_face_at_boundary[{
1582 current_polytope_index,
1584 [current_polytope_index]}] = {true,
1585 neighboring_cell};
1586
1587 const unsigned int n_face =
1589 [current_polytope_index];
1590
1591 is_face_at_boundary[n_face] = true;
1592
1594 [current_polytope_index];
1595
1596 visited_polygonal_neighbors.insert(
1597 std::numeric_limits<unsigned int>::max());
1598 }
1599
1600
1601
1602 if (handler.polytope_cache.visited_cell_and_faces.find(
1603 {cell_index, f}) ==
1604 std::end(handler.polytope_cache.visited_cell_and_faces))
1605 {
1606 handler.polytope_cache
1607 .interface[{
1608 current_polytope_id, current_polytope_id}]
1609 .emplace_back(cell, f);
1610
1611 handler.polytope_cache.visited_cell_and_faces.insert(
1612 {cell_index, f});
1613 }
1614 }
1615 } // loop over faces
1616 } // loop over all cells of agglomerate
1617
1618
1619
1620 if (ghost_counter > 0)
1621 {
1622 const auto parallel_triangulation = dynamic_cast<
1623 const dealii::parallel::TriangulationBase<dim, spacedim> *>(
1624 &(*handler.tria));
1625
1626 const unsigned int n_faces_current_poly =
1627 handler.number_of_agglomerated_faces[current_polytope_index];
1628
1629 // Communicate to neighboring ranks that current_polytope_id has
1630 // a number of faces equal to n_faces_current_poly faces:
1631 // current_polytope_id -> n_faces_current_poly
1632 for (const unsigned int neigh_rank :
1633 parallel_triangulation->ghost_owners())
1634 {
1635 handler.local_n_faces[neigh_rank].emplace(current_polytope_id,
1636 n_faces_current_poly);
1637
1638 handler.local_bdary_info[neigh_rank].emplace(
1639 current_polytope_id, is_face_at_boundary);
1640
1641 handler.local_ghosted_master_id[neigh_rank].emplace(
1642 current_polytope_id, face_to_neigh_id);
1643 }
1644 }
1645 }
1646 };
1647
1648
1649
1650 } // namespace internal
1651} // namespace dealii
1652
1653
1654
1655template class AgglomerationHandler<1>;
1656template void
1658 DynamicSparsityPattern &sparsity_pattern,
1659 const AffineConstraints<double> &constraints,
1660 const bool keep_constrained_dofs,
1661 const types::subdomain_id subdomain_id);
1662
1663template void
1665 TrilinosWrappers::SparsityPattern &sparsity_pattern,
1666 const AffineConstraints<double> &constraints,
1667 const bool keep_constrained_dofs,
1668 const types::subdomain_id subdomain_id);
1669
1670template class AgglomerationHandler<2>;
1671template void
1673 DynamicSparsityPattern &sparsity_pattern,
1674 const AffineConstraints<double> &constraints,
1675 const bool keep_constrained_dofs,
1676 const types::subdomain_id subdomain_id);
1677
1678template void
1680 TrilinosWrappers::SparsityPattern &sparsity_pattern,
1681 const AffineConstraints<double> &constraints,
1682 const bool keep_constrained_dofs,
1683 const types::subdomain_id subdomain_id);
1684
1685template class AgglomerationHandler<3>;
1686template void
1688 DynamicSparsityPattern &sparsity_pattern,
1689 const AffineConstraints<double> &constraints,
1690 const bool keep_constrained_dofs,
1691 const types::subdomain_id subdomain_id);
1692
1693template void
1695 TrilinosWrappers::SparsityPattern &sparsity_pattern,
1696 const AffineConstraints<double> &constraints,
1697 const bool keep_constrained_dofs,
1698 const types::subdomain_id subdomain_id);
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
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const unsigned int degree
size_type size() const
unsigned int size() const
AgglomerationContainer get_agglomerate() const
const AgglomerationIterator< dim, spacedim > neighbor(const unsigned int f) const
const FiniteElement< dim, spacedim > & get_fe() const
DoFHandler< dim, spacedim >::active_cell_iterator as_dof_handler_iterator(const DoFHandler< dim, spacedim > &dof_handler) const
types::subdomain_id subdomain_id() const
unsigned int n_dofs_per_cell() const noexcept
AgglomerationHandler()=default
std::map< types::subdomain_id, std::map< CellId, unsigned int > > local_n_faces
std::unique_ptr< GridTools::Cache< dim, spacedim > > cached_tria
LinearAlgebra::distributed::Vector< float > master_slave_relationships
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Tensor< 1, spacedim > > > > local_normals
const FEValues< dim, spacedim > & reinit(const AgglomerationIterator< dim, spacedim > &polytope) const
std::map< types::subdomain_id, std::map< CellId, CellId > > recv_cell_ids_neigh_cell
FE_Nothing< dim, spacedim > dummy_fe
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv_bdary
std::unique_ptr< hp::FEValues< dim, spacedim > > hp_no_values
bool is_slave_cell(const CellIterator &cell) const
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, bool > > > local_bdary_info
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv
std::unique_ptr< ScratchData > standard_scratch
Quadrature< dim > agglomeration_quad
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< double > > > recv_jxws
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Tensor< 1, spacedim > > > > recv_normals
Quadrature< dim - 1 > agglomeration_face_quad
std::unique_ptr< FEValues< dim, spacedim > > no_values
std::map< types::subdomain_id, std::map< CellId, unsigned int > > recv_n_faces
std::vector< BoundingBox< spacedim > > bboxes
std::unique_ptr< hp::FEFaceValues< dim, spacedim > > hp_no_face_values
std::map< types::subdomain_id, std::map< CellId, BoundingBox< dim > > > recv_ghosted_bbox
bool are_cells_agglomerated(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const typename Triangulation< dim, spacedim >::active_cell_iterator &other_cell) const
Quadrature< dim > agglomerated_quadrature(const AgglomerationContainer &cells, const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell) const
AgglomerationIterator< dim, spacedim > agglomeration_iterator
hp::QCollection< dim - 1 > agglomeration_face_quad_collection
std::map< const typename Triangulation< dim, spacedim >::active_cell_iterator, std::vector< typename Triangulation< dim >::active_face_iterator > > polygon_boundary
ObserverPointer< const Triangulation< dim, spacedim > > tria
std::unordered_map< types::global_cell_index, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > > master2slaves
void create_bounding_box(const AgglomerationContainer &polytope)
std::map< types::global_cell_index, types::global_cell_index > master2polygon
typename AgglomerationIterator< dim, spacedim >::AgglomerationContainer AgglomerationContainer
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Point< spacedim > > > > recv_qpoints
std::pair< const FEValuesBase< dim, spacedim > &, const FEValuesBase< dim, spacedim > & > reinit_interface(const AgglomerationIterator< dim, spacedim > &polytope_in, const AgglomerationIterator< dim, spacedim > &neigh_polytope, const unsigned int local_in, const unsigned int local_outside) const
std::map< types::subdomain_id, std::map< CellId, std::vector< types::global_dof_index > > > recv_ghost_dofs
const Triangulation< dim, spacedim > & get_triangulation() const
void create_agglomeration_sparsity_pattern(SparsityPatternType &sparsity_pattern, const AffineConstraints< Number > &constraints=AffineConstraints< Number >(), const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
const Mapping< dim > & get_mapping() const
std::map< types::global_cell_index, typename Triangulation< dim, spacedim >::active_cell_iterator > master_slave_relationships_iterators
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > get_agglomerate(const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell) const
unsigned int n_agglomerated_faces_per_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
const UpdateFlags internal_agglomeration_face_flags
void distribute_agglomerated_dofs(const FiniteElement< dim > &fe_space)
hp::FECollection< dim, spacedim > dummy_fe_collection
std::vector< types::global_cell_index > number_of_agglomerated_faces
const UpdateFlags internal_agglomeration_flags
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, bool > > > recv_bdary_info
std::unique_ptr< FEFaceValues< dim, spacedim > > no_face_values
const FEValuesBase< dim, spacedim > & reinit_master(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_number, std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > &agglo_isv_ptr) const
std::map< types::subdomain_id, std::map< CellId, std::vector< types::global_dof_index > > > local_ghost_dofs
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > master_cells_container
std::unique_ptr< hp::FECollection< dim, spacedim > > hp_fe_collection
hp::QCollection< dim > agglomeration_quad_collection
bool is_master_cell(const CellIterator &cell) const
IteratorRange< agglomeration_iterator > polytope_iterators() const
void initialize_agglomeration_data(const std::unique_ptr< GridTools::Cache< dim, spacedim > > &cache_tria)
agglomeration_iterator define_agglomerate(const AgglomerationContainer &cells)
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< double > > > local_jxws
void initialize_fe_values(const Quadrature< dim > &cell_quadrature=QGauss< dim >(1), const UpdateFlags &flags=UpdateFlags::update_default, const Quadrature< dim - 1 > &face_quadrature=QGauss< dim - 1 >(1), const UpdateFlags &face_flags=UpdateFlags::update_default)
std::map< types::subdomain_id, std::map< std::pair< CellId, unsigned int >, std::vector< Point< spacedim > > > > local_qpoints
std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > agglomerated_isv_neigh
void define_agglomerate_with_check(const AgglomerationContainer &cells)
std::unique_ptr< ScratchData > agglomerated_scratch
hp::FECollection< dim, spacedim > fe_collection
std::map< types::subdomain_id, std::map< CellId, BoundingBox< dim > > > local_ghosted_bbox
DoFHandler< dim, spacedim > agglo_dh
std::unique_ptr< MappingBox< dim > > box_mapping
std::unique_ptr< FiniteElement< dim > > fe
ObserverPointer< const Mapping< dim, spacedim > > mapping
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, CellId > > > recv_ghosted_master_id
std::map< types::subdomain_id, std::map< CellId, CellId > > local_cell_ids_neigh_cell
internal::PolytopeCache< dim, spacedim > polytope_cache
std::map< types::subdomain_id, std::map< CellId, std::map< unsigned int, CellId > > > local_ghosted_master_id
hp::MappingCollection< dim > mapping_collection
static void setup_master_neighbor_connectivity(const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell, const AgglomerationHandler< dim, spacedim > &handler)
static const FEValuesBase< dim, spacedim > & reinit_master(const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, const unsigned int face_index, std::unique_ptr< NonMatching::FEImmersedSurfaceValues< spacedim > > &agglo_isv_ptr, const AgglomerationHandler< dim, spacedim > &handler)
unsigned int size() const
unsigned int cell_index
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
void make_sparsity_pattern(const DoFHandler< dim, spacedim > &dof_handler, SparsityPatternBase &sparsity_pattern, const AffineConstraints< number > &constraints={}, const bool keep_constrained_dofs=true, const types::subdomain_id subdomain_id=numbers::invalid_subdomain_id)
UpdateFlags
update_normal_vectors
update_JxW_values
update_quadrature_points
IndexSet extract_locally_relevant_dofs(const DoFHandler< dim, spacedim > &dof_handler)
Tensor< 2, dim, Number > w(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
std::map< unsigned int, T > some_to_some(const MPI_Comm comm, const std::map< unsigned int, T > &objects_to_send)
bool job_supports_mpi()
void create_graph_from_agglomerate(const std::vector< typename Triangulation< dim >::active_cell_iterator > &cells, Graph &g)
Definition utils.h:1097
void compute_connected_components(Graph &g, std::vector< std::vector< types::global_cell_index > > &connected_components)
Definition utils.h:1147
constexpr types::subdomain_id invalid_subdomain_id
unsigned int subdomain_id
unsigned int global_cell_index