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
tria_description.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2019 - 2021 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_grid_construction_utilities_h
17#define dealii_grid_construction_utilities_h
18
19#include <deal.II/base/config.h>
20
21#include <deal.II/base/mpi.h>
22
25#include <deal.II/grid/tria.h>
26
28
29
31
32/*------------------------------------------------------------------------*/
33
69template <int structdim>
71{
91 std::vector<unsigned int> vertices;
92
100 union
101 {
111
122 };
123
133
151 CellData(const unsigned int n_vertices =
152 ReferenceCells::get_hypercube<structdim>().n_vertices());
153
157 bool
158 operator==(const CellData<structdim> &other) const;
159
165 template <class Archive>
166 void
167 serialize(Archive &ar, const unsigned int version);
168
169 static_assert(structdim > 0,
170 "The class CellData can only be used for structdim>0.");
171};
172
173
174
230{
238 std::vector<CellData<1>> boundary_lines;
239
254 std::vector<CellData<2>> boundary_quads;
255
264 bool
265 check_consistency(const unsigned int dim) const;
266};
267
268
269template <int structdim>
270template <class Archive>
271void
272CellData<structdim>::serialize(Archive &ar, const unsigned int /*version*/)
273{
274 ar &vertices;
275 ar &material_id;
276 ar &boundary_id;
277 ar &manifold_id;
278}
279
285{
291 {
302 };
303
316 template <int dim>
317 struct CellData
318 {
324 template <class Archive>
325 void
326 serialize(Archive &ar, const unsigned int /*version*/);
327
331 bool
332 operator==(const CellData<dim> &other) const;
333
338
343
348
353
359 std::array<types::manifold_id, GeometryInfo<dim>::lines_per_cell>
361
367 std::array<types::manifold_id,
368 dim == 1 ? 1 : GeometryInfo<3>::quads_per_cell>
370
375 std::vector<std::pair<unsigned int, types::boundary_id>> boundary_ids;
376 };
377
381 template <int dim, int spacedim>
383 {
389 template <class Archive>
390 void
391 serialize(Archive &ar, const unsigned int /*version*/);
392
396 bool
397 operator==(const Description<dim, spacedim> &other) const;
398
402 std::vector<::CellData<dim>> coarse_cells;
403
407 std::vector<Point<spacedim>> coarse_cell_vertices;
408
414 std::vector<types::coarse_cell_id> coarse_cell_index_to_coarse_cell_id;
415
421 std::vector<std::vector<CellData<dim>>> cell_infos;
422
433
438
443 };
444
445
451 namespace Utilities
452 {
518 template <int dim, int spacedim = dim>
521 const ::Triangulation<dim, spacedim> &tria,
522 const MPI_Comm & comm,
525 const unsigned int my_rank_in = numbers::invalid_unsigned_int);
526
547 template <int dim, int spacedim>
554
559 template <int dim, int spacedim>
565 & mg_partitions,
568
626 template <int dim, int spacedim = dim>
629 const std::function<void(::Triangulation<dim, spacedim> &)>
630 & serial_grid_generator,
631 const std::function<void(::Triangulation<dim, spacedim> &,
632 const MPI_Comm &,
633 const unsigned int)> &serial_grid_partitioner,
634 const MPI_Comm & comm,
635 const int group_size = 1,
636 const typename Triangulation<dim, spacedim>::MeshSmoothing smoothing =
640
641 } // namespace Utilities
642
643
644
645 template <int dim>
646 template <class Archive>
647 void
648 CellData<dim>::serialize(Archive &ar, const unsigned int /*version*/)
649 {
650 ar &id;
651 ar &subdomain_id;
652 ar &level_subdomain_id;
653 ar &manifold_id;
654 if (dim >= 2)
655 ar &manifold_line_ids;
656 if (dim >= 3)
657 ar &manifold_quad_ids;
658 ar &boundary_ids;
659 }
660
661
662 template <int dim, int spacedim>
663 template <class Archive>
664 void
666 const unsigned int /*version*/)
667 {
668 ar &coarse_cells;
671 ar &cell_infos;
672 ar &settings;
673 ar &smoothing;
674 }
675
676
677
678 template <int dim>
679 bool
681 {
682 if (this->id != other.id)
683 return false;
684 if (this->subdomain_id != other.subdomain_id)
685 return false;
686 if (this->level_subdomain_id != other.level_subdomain_id)
687 return false;
688 if (this->manifold_id != other.manifold_id)
689 return false;
690 if (dim >= 2 && this->manifold_line_ids != other.manifold_line_ids)
691 return false;
692 if (dim >= 3 && this->manifold_quad_ids != other.manifold_quad_ids)
693 return false;
694 if (this->boundary_ids != other.boundary_ids)
695 return false;
696
697 return true;
698 }
699
700
701
702 template <int dim, int spacedim>
703 bool
705 const Description<dim, spacedim> &other) const
706 {
707 if (this->coarse_cells != other.coarse_cells)
708 return false;
710 return false;
713 return false;
714 if (this->cell_infos != other.cell_infos)
715 return false;
716 if (this->settings != other.settings)
717 return false;
718 if (this->smoothing != other.smoothing)
719 return false;
720
721 return true;
722 }
723} // namespace TriangulationDescription
724
725
727
728#endif
std::array< unsigned int, 4 > binary_type
Definition: cell_id.h:80
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
Point< 3 > vertices[4]
Description< dim, spacedim > create_description_from_triangulation_in_groups(const std::function< void(::Triangulation< dim, spacedim > &)> &serial_grid_generator, const std::function< void(::Triangulation< dim, spacedim > &, const MPI_Comm &, const unsigned int)> &serial_grid_partitioner, const MPI_Comm &comm, const int group_size=1, const typename Triangulation< dim, spacedim >::MeshSmoothing smoothing=::Triangulation< dim, spacedim >::none, const TriangulationDescription::Settings setting=TriangulationDescription::Settings::default_setting)
Description< dim, spacedim > create_description_from_triangulation(const ::Triangulation< dim, spacedim > &tria, const MPI_Comm &comm, const TriangulationDescription::Settings settings=TriangulationDescription::Settings::default_setting, const unsigned int my_rank_in=numbers::invalid_unsigned_int)
static const unsigned int invalid_unsigned_int
Definition: types.h:201
unsigned int manifold_id
Definition: types.h:141
std::vector< unsigned int > vertices
bool operator==(const CellData< structdim > &other) const
types::manifold_id manifold_id
types::material_id material_id
types::boundary_id boundary_id
void serialize(Archive &ar, const unsigned int version)
std::vector< CellData< 2 > > boundary_quads
bool check_consistency(const unsigned int dim) const
std::vector< CellData< 1 > > boundary_lines
std::array< types::manifold_id, dim==1 ? 1 :GeometryInfo< 3 >::quads_per_cell > manifold_quad_ids
std::array< types::manifold_id, GeometryInfo< dim >::lines_per_cell > manifold_line_ids
void serialize(Archive &ar, const unsigned int)
std::vector< std::pair< unsigned int, types::boundary_id > > boundary_ids
bool operator==(const CellData< dim > &other) const
std::vector< std::vector< CellData< dim > > > cell_infos
Triangulation< dim, spacedim >::MeshSmoothing smoothing
void serialize(Archive &ar, const unsigned int)
std::vector<::CellData< dim > > coarse_cells
bool operator==(const Description< dim, spacedim > &other) const
std::vector< Point< spacedim > > coarse_cell_vertices
std::vector< types::coarse_cell_id > coarse_cell_index_to_coarse_cell_id
std::vector< std::pair< unsigned int, Point< spacedim > > > coarse_cell_vertices
const TriangulationDescription::Settings settings
std::vector< types::coarse_cell_id > coarse_cell_index_to_coarse_cell_id
const MPI_Comm & comm
std::vector< std::vector< CellData< dim > > > cell_infos
const ::Triangulation< dim, spacedim > & tria
std::vector<::CellData< dim > > coarse_cells