Reference documentation for deal.II version 9.5.0
\(\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
tria_base.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2008 - 2023 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
16#ifndef dealii_distributed_tria_base_h
17#define dealii_distributed_tria_base_h
18
19
20#include <deal.II/base/config.h>
21
27
28#include <deal.II/grid/tria.h>
29
30#include <functional>
31#include <list>
32#include <set>
33#include <utility>
34#include <vector>
35
36
38
39namespace parallel
40{
78 template <int dim, int spacedim = dim>
80 class TriangulationBase : public ::Triangulation<dim, spacedim>
81 {
82 public:
87 const MPI_Comm mpi_communicator,
88 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
90 const bool check_for_distorted_cells = false);
91
95 virtual ~TriangulationBase() override;
96
100 virtual MPI_Comm
101 get_communicator() const override;
102
106 virtual bool
108
117 virtual void
118 copy_triangulation(
119 const ::Triangulation<dim, spacedim> &old_tria) override;
120
139 unsigned int
140 n_locally_owned_active_cells() const;
141
148 n_global_active_cells() const override;
149
153 virtual std::size_t
154 memory_consumption() const override;
155
156
164 virtual unsigned int
165 n_global_levels() const override;
166
174 locally_owned_subdomain() const override;
175
185 const std::set<types::subdomain_id> &
186 ghost_owners() const;
187
203 const std::set<types::subdomain_id> &
204 level_ghost_owners() const;
205
206 const std::weak_ptr<const Utilities::MPI::Partitioner>
207 global_active_cell_index_partitioner() const override;
208
209 const std::weak_ptr<const Utilities::MPI::Partitioner>
210 global_level_cell_index_partitioner(
211 const unsigned int level) const override;
212
219 virtual std::vector<types::boundary_id>
220 get_boundary_ids() const override;
221
228 virtual std::vector<types::manifold_id>
229 get_manifold_ids() const override;
230
283 void
284 communicate_locally_moved_vertices(
285 const std::vector<bool> &vertex_locally_moved);
286
288 n_global_coarse_cells() const override;
289
290 protected:
297
303
308
314 {
319
325
330
336 unsigned int n_global_levels;
337
342 std::set<types::subdomain_id> ghost_owners;
343
348 std::set<types::subdomain_id> level_ghost_owners;
349
353 std::shared_ptr<const Utilities::MPI::Partitioner>
355
359 std::vector<std::shared_ptr<const Utilities::MPI::Partitioner>>
361
362 NumberCache();
363 };
364
366
370 virtual void
371 update_number_cache();
372
376 void
377 update_reference_cells() override;
378
382 void
383 reset_global_cell_indices();
384 };
385
386
387
435 template <int dim, int spacedim = dim>
438 : public ::parallel::TriangulationBase<dim, spacedim>
439 {
440 public:
445 const MPI_Comm mpi_communicator,
446 const typename ::Triangulation<dim, spacedim>::MeshSmoothing
448 const bool check_for_distorted_cells = false);
449
456 virtual void
457 clear() override;
458
460 typename ::Triangulation<dim, spacedim>::cell_iterator;
461
463 typename ::Triangulation<dim, spacedim>::CellStatus;
464
480 virtual bool
481 has_hanging_nodes() const override;
482
490 virtual void
491 save(const std::string &filename) const = 0;
492
498 virtual void
499 load(const std::string &filename) = 0;
500
507 virtual void
508 load(const std::string &filename, const bool autopartition) = 0;
509
618 unsigned int
619 register_data_attach(
620 const std::function<std::vector<char>(const cell_iterator &,
621 const CellStatus)> &pack_callback,
622 const bool returns_variable_size_data);
623
672 void
673 notify_ready_to_unpack(
674 const unsigned int handle,
675 const std::function<
676 void(const cell_iterator &,
677 const CellStatus,
678 const boost::iterator_range<std::vector<char>::const_iterator> &)>
679 &unpack_callback);
680
681 protected:
689 void
690 save_attached_data(const unsigned int global_first_cell,
691 const unsigned int global_num_cells,
692 const std::string &filename) const;
693
702 void
703 load_attached_data(const unsigned int global_first_cell,
704 const unsigned int global_num_cells,
705 const unsigned int local_num_cells,
706 const std::string &filename,
707 const unsigned int n_attached_deserialize_fixed,
708 const unsigned int n_attached_deserialize_variable);
709
720 virtual void
722
728 using cell_relation_t = typename std::pair<cell_iterator, CellStatus>;
729
735 std::vector<cell_relation_t> local_cell_relations;
736
743 {
749
755
756 using pack_callback_t = std::function<std::vector<char>(
757 typename ::Triangulation<dim, spacedim>::cell_iterator,
758 typename ::Triangulation<dim, spacedim>::CellStatus)>;
759
764 std::vector<pack_callback_t> pack_callbacks_fixed;
765 std::vector<pack_callback_t> pack_callbacks_variable;
766 };
767
769
778 {
779 public:
780 DataTransfer(const MPI_Comm mpi_communicator);
781
791 void
792 pack_data(const std::vector<cell_relation_t> &cell_relations,
793 const std::vector<typename CellAttachedData::pack_callback_t>
794 &pack_callbacks_fixed,
795 const std::vector<typename CellAttachedData::pack_callback_t>
796 &pack_callbacks_variable);
797
798
799
807 void
808 unpack_cell_status(std::vector<cell_relation_t> &cell_relations) const;
809
822 void
823 unpack_data(
824 const std::vector<cell_relation_t> &cell_relations,
825 const unsigned int handle,
826 const std::function<void(
827 const typename ::Triangulation<dim, spacedim>::cell_iterator &,
828 const typename ::Triangulation<dim, spacedim>::CellStatus &,
829 const boost::iterator_range<std::vector<char>::const_iterator> &)>
830 &unpack_callback) const;
831
846 void
847 save(const unsigned int global_first_cell,
848 const unsigned int global_num_cells,
849 const std::string &filename) const;
850
869 void
870 load(const unsigned int global_first_cell,
871 const unsigned int global_num_cells,
872 const unsigned int local_num_cells,
873 const std::string &filename,
874 const unsigned int n_attached_deserialize_fixed,
875 const unsigned int n_attached_deserialize_variable);
876
883 void
884 clear();
885
890
901 std::vector<unsigned int> sizes_fixed_cumulative;
902
907 std::vector<char> src_data_fixed;
908 std::vector<char> dest_data_fixed;
909
914 std::vector<int> src_sizes_variable;
915 std::vector<int> dest_sizes_variable;
916 std::vector<char> src_data_variable;
917 std::vector<char> dest_data_variable;
918
919 private:
921 };
922
924 };
925
926} // namespace parallel
927
929
930#endif
std::vector< unsigned int > sizes_fixed_cumulative
Definition tria_base.h:901
virtual void load(const std::string &filename)=0
virtual void load(const std::string &filename, const bool autopartition)=0
typename ::Triangulation< dim, spacedim >::CellStatus CellStatus
Definition tria_base.h:463
std::vector< cell_relation_t > local_cell_relations
Definition tria_base.h:735
virtual void save(const std::string &filename) const =0
typename ::Triangulation< dim, spacedim >::cell_iterator cell_iterator
Definition tria_base.h:460
typename std::pair< cell_iterator, CellStatus > cell_relation_t
Definition tria_base.h:728
types::subdomain_id n_subdomains
Definition tria_base.h:307
const MPI_Comm mpi_communicator
Definition tria_base.h:296
virtual bool is_multilevel_hierarchy_constructed() const =0
types::subdomain_id my_subdomain
Definition tria_base.h:302
#define DEAL_II_DEPRECATED
Definition config.h:172
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:160
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
unsigned int level
Definition grid_out.cc:4618
std::vector< pack_callback_t > pack_callbacks_fixed
Definition tria_base.h:764
std::vector< pack_callback_t > pack_callbacks_variable
Definition tria_base.h:765
std::function< std::vector< char >(typename ::Triangulation< dim, spacedim >::cell_iterator, typename ::Triangulation< dim, spacedim >::CellStatus)> pack_callback_t
Definition tria_base.h:758
std::set< types::subdomain_id > level_ghost_owners
Definition tria_base.h:348
types::global_cell_index n_global_active_cells
Definition tria_base.h:324
std::set< types::subdomain_id > ghost_owners
Definition tria_base.h:342
std::shared_ptr< const Utilities::MPI::Partitioner > active_cell_index_partitioner
Definition tria_base.h:354
std::vector< std::shared_ptr< const Utilities::MPI::Partitioner > > level_cell_index_partitioners
Definition tria_base.h:360
types::coarse_cell_id number_of_global_coarse_cells
Definition tria_base.h:329