16 #include <deal.II/base/memory_consumption.h> 17 #include <deal.II/base/smartpointer.h> 19 #include <deal.II/distributed/shared_tria.h> 20 #include <deal.II/distributed/tria.h> 22 #include <deal.II/dofs/dof_accessor.h> 23 #include <deal.II/dofs/dof_handler.h> 25 #include <deal.II/fe/fe.h> 27 #include <deal.II/grid/intergrid_map.h> 28 #include <deal.II/grid/tria.h> 29 #include <deal.II/grid/tria_accessor.h> 30 #include <deal.II/grid/tria_iterator.h> 33 DEAL_II_NAMESPACE_OPEN
36 template <
class MeshType>
38 : source_grid(nullptr, typeid(*this).name())
39 , destination_grid(nullptr, typeid(*this).name())
44 template <
class MeshType>
47 const MeshType &destination_grid)
53 this->source_grid = &source_grid;
54 this->destination_grid = &destination_grid;
58 const unsigned int n_levels = source_grid.get_triangulation().n_levels();
59 mapping.resize(n_levels);
60 for (
unsigned int level = 0; level < n_levels; ++level)
69 unsigned int n_cells = 0;
71 endc = source_grid.end(level);
72 for (; cell != endc; ++cell)
73 if (static_cast<unsigned int>(cell->index()) > n_cells)
74 n_cells = cell->index();
79 mapping[level].resize(n_cells + 1, destination_grid.end());
89 dst_cell = destination_grid.begin(0), endc = source_grid.end(0);
90 for (; src_cell != endc; ++src_cell, ++dst_cell)
91 set_mapping(src_cell, dst_cell);
95 Assert(dst_cell == destination_grid.end(0), ExcIncompatibleGrids());
100 template <
class MeshType>
106 mapping[src_cell->level()][src_cell->index()] = dst_cell;
110 if (src_cell->has_children() && dst_cell->has_children())
112 Assert(src_cell->n_children() ==
115 Assert(dst_cell->n_children() ==
118 Assert(src_cell->refinement_case() == dst_cell->refinement_case(),
120 for (
unsigned int c = 0;
121 c < GeometryInfo<MeshType::dimension>::max_children_per_cell;
123 set_mapping(src_cell->child(c), dst_cell->child(c));
125 else if (src_cell->has_children() && !dst_cell->has_children())
130 for (
unsigned int c = 0; c < src_cell->n_children(); ++c)
131 set_entries_to_cell(src_cell->child(c), dst_cell);
139 template <
class MeshType>
145 mapping[src_cell->level()][src_cell->index()] = dst_cell;
149 if (src_cell->has_children())
150 for (
unsigned int c = 0; c < src_cell->n_children(); ++c)
151 set_entries_to_cell(src_cell->child(c), dst_cell);
155 template <
class MeshType>
160 ExcInvalidKey(source_cell));
161 Assert(source_cell->level() <=
static_cast<int>(mapping.size()),
162 ExcInvalidKey(source_cell));
163 Assert(source_cell->index() <=
164 static_cast<int>(mapping[source_cell->level()].size()),
165 ExcInvalidKey(source_cell));
167 return mapping[source_cell->level()][source_cell->index()];
172 template <
class MeshType>
177 source_grid =
nullptr;
178 destination_grid =
nullptr;
183 template <
class MeshType>
192 template <
class MeshType>
196 return *destination_grid;
201 template <
class MeshType>
213 #include "intergrid_map.inst" 215 DEAL_II_NAMESPACE_CLOSE
cell_iterator operator[](const cell_iterator &source_cell) const
const MeshType & get_destination_grid() const
void set_entries_to_cell(const cell_iterator &src_cell, const cell_iterator &dst_cell)
typename MeshType::cell_iterator cell_iterator
const MeshType & get_source_grid() const
std::size_t memory_consumption() const
#define Assert(cond, exc)
void make_mapping(const MeshType &source_grid, const MeshType &destination_grid)
static ::ExceptionBase & ExcNotImplemented()
Iterator points to a valid object.
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)
void set_mapping(const cell_iterator &src_cell, const cell_iterator &dst_cell)