PolyDEAL
 
Loading...
Searching...
No Matches
agglomeration_accessor.h
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#ifndef agglomeration_accessor_h
15#define agglomeration_accessor_h
16
17#include <deal.II/base/config.h>
18
21
23
24#include <vector>
25
26using namespace dealii;
27
28
29// Forward declarations
30#ifndef DOXYGEN
31template <int, int>
33template <int, int>
35#endif
36
37
41template <int dim, int spacedim = dim>
43{
44public:
49 std::vector<typename Triangulation<dim, spacedim>::active_cell_iterator>;
50
51
55 void
56 get_dof_indices(std::vector<types::global_dof_index> &) const;
57
64 unsigned int
65 n_faces() const;
66
70 unsigned int
72
77 neighbor(const unsigned int f) const;
78
83 unsigned int
84 neighbor_of_agglomerated_neighbor(const unsigned int f) const;
85
95 bool
96 at_boundary(const unsigned int f) const;
97
101 const std::vector<typename Triangulation<dim>::active_face_iterator> &
103
108 double
109 volume() const;
110
114 double
115 diameter() const;
116
122
128 const BoundingBox<dim> &
130
135 index() const;
136
145
150 unsigned int
152
153 /* Returns true if this polygon is owned by the current processor. On a serial
154 * Triangulation this returs always true, but may yield false for a
155 * parallel::distributed::Triangulation.
156 */
157 bool
159
165 CellId
166 id() const;
167
174
178 inline const std::vector<types::global_cell_index> &
179 children() const;
180
181private:
187
197
203 const CellId &cell_id,
205
210
211
216
221
226
231
236
241 bool
243
247 bool
249
253 void
255
259 void
261
266 get_slaves() const;
267
268 unsigned int
271 const;
272
273 template <int, int>
275};
276
277
278
279template <int dim, int spacedim>
280unsigned int
283{
284 unsigned int n_neighbors = 0;
285 for (const auto &f : cell->face_indices())
286 {
287 const auto &neighboring_cell = cell->neighbor(f);
288 if ((cell->face(f)->at_boundary()) ||
289 (neighboring_cell->is_active() &&
290 !handler->are_cells_agglomerated(cell, neighboring_cell)))
291 {
292 ++n_neighbors;
293 }
294 }
295 return n_neighbors;
296}
297
298
299
300template <int dim, int spacedim>
301unsigned int
303{
304 Assert(!handler->is_slave_cell(master_cell),
305 ExcMessage("You cannot pass a slave cell."));
306 return handler->number_of_agglomerated_faces[present_index];
307}
308
309
310
311template <int dim, int spacedim>
314{
315 if (!at_boundary(f))
316 {
317 if (master_cell->is_ghost())
318 {
319 // The following path is needed when the present function is called
320 // from neighbor_of_neighbor()
321
322 const unsigned int sender_rank = master_cell->subdomain_id();
323
324 const CellId &master_id_ghosted_neighbor =
325 handler->recv_ghosted_master_id.at(sender_rank)
326 .at(present_id)
327 .at(f);
328
329 // Use the id of the master cell to uniquely identify the neighboring
330 // agglomerate
331
332 return {master_cell,
333 master_id_ghosted_neighbor,
334 handler}; // dummy master?
335 }
336
337 const types::global_cell_index polytope_index =
338 handler->master2polygon.at(master_cell->active_cell_index());
339
340 // #ifdef AGGLO_DEBUG
341 // std::cout << "Poly index: " << polytope_index << std::endl;
342 // std::cout << "Face index " << f << std::endl;
343
344 // if (Utilities::MPI::this_mpi_process(handler->communicator) == 1)
345 // for (const auto &[key, value] :
346 // handler->polytope_cache.cell_face_at_boundary)
347 // {
348 // std::cout << "p.f = " << key.first << std::endl;
349 // std::cout << "p.s = " << key.second << std::endl;
350 // }
351 // #endif
352
353 const auto &neigh =
354 handler->polytope_cache.cell_face_at_boundary.at({polytope_index, f})
355 .second;
356
357
358 if (neigh->is_locally_owned())
359 {
361 *neigh, &(handler->agglo_dh));
362 return {cell_dh, handler};
363 }
364 else
365 {
366 // #ifdef AGGLO_DEBUG
367 // std::cout << "before calling ghosted_masteer_id" <<
368 // std::endl; std::cout << "present_id = " << present_id <<
369 // std::endl; std::cout << "che master cell? " <<
370 // master_cell << std::endl;
371
372 // if
373 // (Utilities::MPI::this_mpi_process(handler->communicator)
374 // == 1)
375 // for (const auto &[key, value] :
376 // handler->polytope_cache.ghosted_master_id)
377 // {
378 // std::cout << "poly_id = " << key.first <<
379 // std::endl; std::cout << "face_idx = " << key.second
380 // << std::endl;
381 // }
382
383 // #endif
384 // Get master_id from the neighboring ghost polytope. This uniquely
385 // identifies the neighboring polytope among all processors.
386 const CellId &master_id_neighbor =
387 handler->polytope_cache.ghosted_master_id.at({present_id, f});
388
389 // Use the id of the master cell to uniquely identify the neighboring
390 // agglomerate
391 return {neigh, master_id_neighbor, handler};
392 }
393 }
394 else
395 {
396 return {};
397 }
398}
399
400
401
402template <int dim, int spacedim>
403unsigned int
405 const unsigned int f) const
406{
407 // First, make sure it's not a boundary face.
408 if (!at_boundary(f))
409 {
410 const auto &neigh_polytope =
411 neighbor(f); // returns the neighboring master and id
412
413 AssertThrow(neigh_polytope.state() == IteratorState::valid,
414 ExcInternalError());
415
416 unsigned int n_faces_agglomerated_neighbor;
417
418 // if it is locally owned, retrieve the number of faces
419 if (neigh_polytope->is_locally_owned())
420 {
421 n_faces_agglomerated_neighbor = neigh_polytope->n_faces();
422 }
423 else
424 {
425 // The neighboring polytope is not locally owned. We need to get the
426 // number of its faces from the neighboring rank.
427
428 // First, retrieve the CellId of the neighboring polytope.
429 const CellId &master_id_neighbor = neigh_polytope->id();
430
431 // Then, get the neighboring rank
432 const unsigned int sender_rank = neigh_polytope->subdomain_id();
433
434 // From the neighboring rank, use the CellId of the neighboring
435 // polytope to get the number of its faces.
436 n_faces_agglomerated_neighbor =
437 handler->recv_n_faces.at(sender_rank).at(master_id_neighbor);
438 }
439
440
441 // Loop over all faces of neighboring agglomerate
442 for (unsigned int f_out = 0; f_out < n_faces_agglomerated_neighbor;
443 ++f_out)
444 {
445 // Check if same CellId
446 if (neigh_polytope->neighbor(f_out).state() == IteratorState::valid)
447 if (neigh_polytope->neighbor(f_out)->id() == present_id)
448 return f_out;
449 }
451 }
452 else
453 {
454 // Face is at boundary
456 }
457}
458
459// ------------------------------ inline functions -------------------------
460
461template <int dim, int spacedim>
464
465
466
467template <int dim, int spacedim>
471{
472 handler = const_cast<AgglomerationHandler<dim, spacedim> *>(ah);
473 if (&(*handler->master_cells_container.end()) == std::addressof(cell))
474 {
475 present_index = handler->master_cells_container.size();
476 master_cell = *handler->master_cells_container.end();
477 present_id = CellId(); // invalid id (TODO)
478 present_subdomain_id = numbers::invalid_subdomain_id;
479 }
480 else
481 {
482 present_index = handler->master2polygon.at(cell->active_cell_index());
483 master_cell = cell;
484 present_id = master_cell->id();
485 present_subdomain_id = master_cell->subdomain_id();
486 }
487}
488
489
490
491template <int dim, int spacedim>
494 const CellId &master_cell_id,
496{
497 Assert(neigh_cell->is_ghost(), ExcInternalError());
498 // neigh_cell is ghosted
499
500 handler = const_cast<AgglomerationHandler<dim, spacedim> *>(ah);
501 master_cell = neigh_cell;
502 present_index = numbers::invalid_unsigned_int;
503 // neigh_cell is ghosted, use the CellId of that agglomerate
504 present_id = master_cell_id;
505 present_subdomain_id = master_cell->subdomain_id();
506}
507
508
509
510template <int dim, int spacedim>
511inline void
513 std::vector<types::global_dof_index> &dof_indices) const
514{
515 Assert(dof_indices.size() > 0,
516 ExcMessage(
517 "The vector of DoFs indices must be already properly resized."));
518 if (is_locally_owned())
519 {
520 // Forward the call to the master cell
521 typename DoFHandler<dim, spacedim>::cell_iterator master_cell_dh(
522 *master_cell, &(handler->agglo_dh));
523 master_cell_dh->get_dof_indices(dof_indices);
524 }
525 else
526 {
527 const std::vector<types::global_dof_index> &recv_dof_indices =
528 handler->recv_ghost_dofs.at(present_subdomain_id).at(present_id);
529
530 std::copy(recv_dof_indices.cbegin(),
531 recv_dof_indices.cend(),
532 dof_indices.begin());
533 }
534}
535
536
537
538template <int dim, int spacedim>
541{
542 auto agglomeration = get_slaves();
543 agglomeration.push_back(master_cell);
544 return agglomeration;
545}
546
547
548
549template <int dim, int spacedim>
550inline const std::vector<typename Triangulation<dim>::active_face_iterator> &
552{
553 return handler->polygon_boundary[master_cell];
554}
555
556
557
558template <int dim, int spacedim>
559inline double
561{
562 Assert(!handler->is_slave_cell(master_cell),
563 ExcMessage("The present function cannot be called for slave cells."));
564
565 if (handler->is_master_cell(master_cell))
566 {
567 // Get the bounding box associated with the master cell
568 const auto &bdary_pts =
569 handler->bboxes[present_index].get_boundary_points();
570 return (bdary_pts.second - bdary_pts.first).norm();
571 }
572 else
573 {
574 // Standard deal.II way to get the measure of a cell.
575 return master_cell->diameter();
576 }
577}
578
579
580
581template <int dim, int spacedim>
582inline const BoundingBox<dim> &
584{
585 if (is_locally_owned())
586 return handler->bboxes[present_index];
587 else
588 return handler->recv_ghosted_bbox.at(present_subdomain_id).at(present_id);
589}
590
591
592
593template <int dim, int spacedim>
594inline double
596{
597 Assert(!handler->is_slave_cell(master_cell),
598 ExcMessage("The present function cannot be called for slave cells."));
599
600 if (handler->is_master_cell(master_cell))
601 {
602 return handler->bboxes[present_index].volume();
603 }
604 else
605 {
606 return master_cell->measure();
607 }
608}
609
610
611
612template <int dim, int spacedim>
613inline void
615{
616 // Increment the present index and update the polytope
617 ++present_index;
618
619 // Make sure not to query the CellId if it's past the last
620 if (present_index < handler->master_cells_container.size())
621 {
622 master_cell = handler->master_cells_container[present_index];
623 present_id = master_cell->id();
624 present_subdomain_id = master_cell->subdomain_id();
625 }
626}
627
628
629
630template <int dim, int spacedim>
631inline void
633{
634 // Decrement the present index and update the polytope
635 --present_index;
636 master_cell = handler->master_cells_container[present_index];
637 present_id = master_cell->id();
638}
639
640
641template <int dim, int spacedim>
642inline bool
644 const AgglomerationAccessor<dim, spacedim> &other) const
645{
646 return present_index == other.present_index;
647}
648
649template <int dim, int spacedim>
650inline bool
652 const AgglomerationAccessor<dim, spacedim> &other) const
653{
654 return !(*this == other);
655}
656
657
658
659template <int dim, int spacedim>
662{
663 return present_index;
664}
665
666
667
668template <int dim, int spacedim>
671 const DoFHandler<dim, spacedim> &dof_handler) const
672{
673 // Forward the call to the master cell using the right DoFHandler.
674 return master_cell->as_dof_handler_iterator(dof_handler);
675}
676
677
678
679template <int dim, int spacedim>
680inline const typename AgglomerationAccessor<dim,
681 spacedim>::AgglomerationContainer &
683{
684 return handler->master2slaves.at(master_cell->active_cell_index());
685}
686
687
688
689template <int dim, int spacedim>
690inline unsigned int
692{
693 AssertThrow(get_agglomerate().size() > 0, ExcMessage("Empty agglomeration."));
694 return get_agglomerate().size();
695}
696
697
698
699template <int dim, int spacedim>
700unsigned int
702{
703 const auto &agglomeration = get_agglomerate();
704 unsigned int n_neighbors = 0;
705 for (const auto &cell : agglomeration)
706 n_neighbors += n_agglomerated_faces_per_cell(cell);
707 return n_neighbors;
708}
709
710
711
712template <int dim, int spacedim>
713inline bool
715{
716 if (master_cell->is_ghost())
717 {
718 const unsigned int sender_rank = master_cell->subdomain_id();
719 // #ifdef AGGLO_DEBUG
720 // const auto & map_recv =
721 // handler->recv_bdary_info.at(sender_rank).at(present_id);
722 // std::cout << "size of map " << map_recv.size() << std::endl;
723 // for (const auto &[face_idx, boundary_or_not] : map_recv)
724 // {
725 // // std::cout << "id: " << face_idx << std::endl;
726 // std::cout << "face: " << face_idx << std::endl;
727 // std::cout << "boundary or not: " << boundary_or_not <<
728 // std::endl;
729 // }
730
731 // std::cout << "il present_id e': " << present_id << "and size "
732 // << map_recv.size() << std::endl;
733 // #endif
734
735 return handler->recv_bdary_info.at(sender_rank).at(present_id).at(f);
736 }
737 else
738 {
739 Assert(!handler->is_slave_cell(master_cell),
740 ExcMessage(
741 "This function should not be called for a slave cell."));
742
743
745 *master_cell, &(handler->agglo_dh));
746 return handler->at_boundary(cell_dh, f);
747 }
748}
749
750
751
752template <int dim, int spacedim>
753inline bool
755{
756 return master_cell->is_locally_owned();
757}
758
759
760
761template <int dim, int spacedim>
762inline CellId
764{
765 return present_id;
766}
767
768
769
770template <int dim, int spacedim>
773{
774 return present_subdomain_id;
775}
776
777template <int dim, int spacedim>
778inline const std::vector<types::global_cell_index> &
780{
781 Assert(!handler->parent_child_info.empty(), ExcInternalError());
782 return handler->parent_child_info.at(
783 {present_index, handler->present_extraction_level});
784}
785
786
787#endif
cell_iterator end() const
unsigned int n_agglomerated_faces() const
AgglomerationContainer get_agglomerate() const
types::global_cell_index index() const
~AgglomerationAccessor()=default
std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > AgglomerationContainer
unsigned int n_agglomerated_faces_per_cell(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
unsigned int n_background_cells() const
unsigned int neighbor_of_agglomerated_neighbor(const unsigned int f) const
AgglomerationAccessor(const typename Triangulation< dim, spacedim >::active_cell_iterator &master_cell, const AgglomerationHandler< dim, spacedim > *ah)
const std::vector< typename Triangulation< dim >::active_face_iterator > & polytope_boundary() const
const AgglomerationIterator< dim, spacedim > neighbor(const unsigned int f) const
bool at_boundary(const unsigned int f) const
const BoundingBox< dim > & get_bounding_box() const
const std::vector< types::global_cell_index > & children() const
bool operator!=(const AgglomerationAccessor< dim, spacedim > &other) const
types::global_cell_index present_index
bool operator==(const AgglomerationAccessor< dim, spacedim > &other) const
AgglomerationHandler< dim, spacedim > * handler
DoFHandler< dim, spacedim >::active_cell_iterator as_dof_handler_iterator(const DoFHandler< dim, spacedim > &dof_handler) const
void get_dof_indices(std::vector< types::global_dof_index > &) const
types::subdomain_id subdomain_id() const
AgglomerationAccessor(const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const CellId &cell_id, const AgglomerationHandler< dim, spacedim > *ah)
const AgglomerationContainer & get_slaves() const
types::subdomain_id present_subdomain_id
Triangulation< dim, spacedim >::active_cell_iterator master_cell
unsigned int n_faces() const
const Triangulation< dim, spacedim >::active_cell_iterator & master_cell() const
typename AgglomerationAccessor< dim, spacedim >::AgglomerationContainer AgglomerationContainer
Point< 2 > second
#define Assert(cond, exc)
#define AssertThrow(cond, exc)
const types::subdomain_id invalid_subdomain_id
static const unsigned int invalid_unsigned_int
unsigned int subdomain_id
unsigned int global_cell_index