Reference documentation for deal.II version GIT relicensing-487-ge9eb5ab491 2024-04-25 07:20:02+00:00
\(\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// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2020 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_grid_construction_utilities_h
16#define dealii_grid_construction_utilities_h
17
18#include <deal.II/base/config.h>
19
22
25#include <deal.II/grid/tria.h>
26
27
29#ifndef DOXYGEN
30namespace LinearAlgebra
31{
32 namespace distributed
33 {
34 template <typename Number, typename MemorySpace>
35 class Vector;
36 }
37} // namespace LinearAlgebra
38#endif
39
40/*------------------------------------------------------------------------*/
41
77template <int structdim>
79{
99 std::vector<unsigned int> vertices;
100
108 union
109 {
119
130 };
131
141
159 CellData(const unsigned int n_vertices =
160 ReferenceCells::get_hypercube<structdim>().n_vertices());
161
165 bool
166 operator==(const CellData<structdim> &other) const;
167
173 template <class Archive>
174 void
175 serialize(Archive &ar, const unsigned int version);
176
177 static_assert(structdim > 0,
178 "The class CellData can only be used for structdim>0.");
179};
180
181
182
238{
246 std::vector<CellData<1>> boundary_lines;
247
262 std::vector<CellData<2>> boundary_quads;
263
272 bool
273 check_consistency(const unsigned int dim) const;
274};
275
276
277template <int structdim>
278template <class Archive>
279void
280CellData<structdim>::serialize(Archive &ar, const unsigned int /*version*/)
281{
282 ar &vertices;
283 ar &material_id;
284 ar &boundary_id;
285 ar &manifold_id;
286}
287
303{
321
338 template <int dim>
339 struct CellData
340 {
344 CellData();
345
351 template <class Archive>
352 void
353 serialize(Archive &ar, const unsigned int /*version*/);
354
358 bool
359 operator==(const CellData<dim> &other) const;
360
365
370
375
380
386 std::array<types::manifold_id, GeometryInfo<dim>::lines_per_cell>
388
394 std::array<types::manifold_id,
395 dim == 1 ? 1 : GeometryInfo<3>::quads_per_cell>
397
402 std::vector<std::pair<unsigned int, types::boundary_id>> boundary_ids;
403 };
404
405
412 template <int dim, int spacedim = dim>
414 {
418 Description();
419
425 template <class Archive>
426 void
427 serialize(Archive &ar, const unsigned int /*version*/);
428
432 bool
433 operator==(const Description<dim, spacedim> &other) const;
434
438 std::vector<::CellData<dim>> coarse_cells;
439
443 std::vector<Point<spacedim>> coarse_cell_vertices;
444
450 std::vector<types::coarse_cell_id> coarse_cell_index_to_coarse_cell_id;
451
457 std::vector<std::vector<CellData<dim>>> cell_infos;
458
469
474
479 };
480
481
487 namespace Utilities
488 {
575 template <int dim, int spacedim = dim>
578 const ::Triangulation<dim, spacedim> &tria,
579 const MPI_Comm comm,
582 const unsigned int my_rank_in = numbers::invalid_unsigned_int);
583
613 template <int dim, int spacedim>
618 &partition,
621
629 template <int dim, int spacedim>
634 &partition,
635 const std::vector<
637 &mg_partitions,
640
698 template <int dim, int spacedim = dim>
701 const std::function<void(::Triangulation<dim, spacedim> &)>
702 &serial_grid_generator,
703 const std::function<void(::Triangulation<dim, spacedim> &,
704 const MPI_Comm,
705 const unsigned int)> &serial_grid_partitioner,
706 const MPI_Comm comm,
707 const int group_size = 1,
708 const typename Triangulation<dim, spacedim>::MeshSmoothing smoothing =
712
713 } // namespace Utilities
714
715
716
717 template <int dim>
719 : subdomain_id(numbers::invalid_subdomain_id)
720 , level_subdomain_id(numbers::invalid_subdomain_id)
721 , manifold_id(numbers::flat_manifold_id)
722 {
723 std::fill(id.begin(), id.end(), numbers::invalid_unsigned_int);
724 std::fill(manifold_line_ids.begin(),
725 manifold_line_ids.end(),
727 std::fill(manifold_quad_ids.begin(),
728 manifold_quad_ids.end(),
730 }
731
732
733
734 template <int dim>
735 template <class Archive>
736 void
737 CellData<dim>::serialize(Archive &ar, const unsigned int /*version*/)
738 {
739 ar &id;
740 ar &subdomain_id;
741 ar &level_subdomain_id;
742 ar &manifold_id;
743 if (dim >= 2)
744 ar &manifold_line_ids;
745 if (dim >= 3)
746 ar &manifold_quad_ids;
747 ar &boundary_ids;
748 }
749
750
751
752 template <int dim, int spacedim>
754 : comm(MPI_COMM_NULL)
755 , settings(Settings::default_setting)
756 , smoothing(Triangulation<dim, spacedim>::MeshSmoothing::none)
757 {}
758
759
760
761 template <int dim, int spacedim>
762 template <class Archive>
763 void
765 const unsigned int /*version*/)
766 {
767 ar &coarse_cells;
770 ar &cell_infos;
771 ar &settings;
772 ar &smoothing;
773 }
774
775
776
777 template <int dim>
778 bool
780 {
781 if (this->id != other.id)
782 return false;
783 if (this->subdomain_id != other.subdomain_id)
784 return false;
785 if (this->level_subdomain_id != other.level_subdomain_id)
786 return false;
787 if (this->manifold_id != other.manifold_id)
788 return false;
789 if (dim >= 2 && this->manifold_line_ids != other.manifold_line_ids)
790 return false;
791 if (dim >= 3 && this->manifold_quad_ids != other.manifold_quad_ids)
792 return false;
793 if (this->boundary_ids != other.boundary_ids)
794 return false;
795
796 return true;
797 }
798
799
800
801 template <int dim, int spacedim>
802 bool
804 const Description<dim, spacedim> &other) const
805 {
806 if (this->coarse_cells != other.coarse_cells)
807 return false;
809 return false;
812 return false;
813 if (this->cell_infos != other.cell_infos)
814 return false;
815 if (this->settings != other.settings)
816 return false;
817 if (this->smoothing != other.smoothing)
818 return false;
819
820 return true;
821 }
822} // namespace TriangulationDescription
823
824
826
827#endif
std::array< unsigned int, 4 > binary_type
Definition cell_id.h:80
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
Point< 3 > vertices[4]
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)
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)
static const unsigned int invalid_unsigned_int
Definition types.h:220
const types::manifold_id flat_manifold_id
Definition types.h:325
unsigned int manifold_id
Definition types.h:156
*braid_SplitCommworld & comm
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
std::vector< types::coarse_cell_id > coarse_cell_index_to_coarse_cell_id
std::vector< std::vector< CellData< dim > > > cell_infos
std::vector<::CellData< dim > > coarse_cells