16#ifndef dealii_grid_generator_h
17#define dealii_grid_generator_h
93 template <
int dim,
int spacedim>
96 const double left = 0.,
97 const double right = 1.,
137 template <
int dim,
int spacedim>
169 template <
int dim,
int spacedim>
172 const unsigned int repetitions,
173 const double left = 0.,
174 const double right = 1.,
212 template <
int dim,
int spacedim>
268 template <
int dim,
int spacedim>
271 const std::vector<unsigned int> &repetitions,
294 const std::vector<std::vector<double>> &step_sizes,
312 const std::vector<std::vector<double>> &spacing,
341 template <
int dim,
int spacedim>
344 const std::vector<unsigned int> &holes);
392 const double inner_radius = 0.4,
393 const double outer_radius = 1.,
394 const double pad_bottom = 2.,
395 const double pad_top = 2.,
396 const double pad_left = 1.,
397 const double pad_right = 1.,
402 const unsigned int n_slices = 2,
494 const double shell_region_width = 0.03,
495 const unsigned int n_shells = 2,
496 const double skewness = 2.0,
519 template <
int dim,
int spacedim>
582 const unsigned int n_subdivisions,
598 const unsigned int (&n_subdivisions)[dim],
600 const unsigned int *n_subdivisions,
628 template <
int dim,
int spacedim>
633 const std::vector<unsigned int> &subdivisions = {},
655 const double left = 0.,
656 const double right = 1.,
657 const double thickness = 1.,
740 const double radius = 1.,
741 const bool attach_spherical_manifold_on_boundary_cells =
false);
779 const double radius = 1.);
802 const bool rotate_left_square,
803 const bool rotate_right_square);
827 const bool face_orientation,
828 const bool face_flip,
829 const bool face_rotation,
830 const bool manipulate_left_cube);
858 template <
int spacedim>
861 const double radius = 1.);
900 const double radius = 1.);
918 const double radius = 1.);
946 const double radius = 1.,
947 const double half_length = 1.);
986 const unsigned int x_subdivisions,
987 const double radius = 1.,
988 const double half_length = 1.);
1024 const double radius_0 = 1.0,
1025 const double radius_1 = 0.5,
1026 const double half_length = 1.0);
1054 template <
int dim,
int spacedim>
1057 const std::vector<unsigned int> &sizes,
1058 const bool colorize_cells =
false);
1098 const double left = -1.,
1099 const double right = 1.,
1133 template <
int dim,
int spacedim>
1136 const std::vector<unsigned int> &repetitions,
1139 const std::vector<int> & n_cells_to_remove);
1163 const double left = 0.,
1164 const double right = 1.,
1228 const double inner_radius,
1229 const double outer_radius,
1230 const unsigned int n_cells = 0,
1271 const double inner_radius,
1272 const double outer_radius,
1308 const double inner_radius,
1309 const double outer_radius,
1310 const unsigned int n_cells = 0,
1343 const double inner_radius,
1344 const double outer_radius,
1345 const unsigned int n_cells = 0,
1377 const double length,
1378 const double inner_radius,
1379 const double outer_radius,
1380 const unsigned int n_radial_cells = 0,
1381 const unsigned int n_axial_cells = 0);
1432 template <
int dim,
int spacedim>
1437 const unsigned int n_cells_toroidal = 6,
1472 const double inner_radius = .25,
1473 const double outer_radius = .5,
1474 const double L = .5,
1475 const unsigned int repetitions = 1,
1552 const double inner_radius = 0.125,
1553 const double outer_radius = 0.25,
1554 const unsigned int n_shells = 1,
1555 const double skewness = 0.1,
1556 const unsigned int n_cells_per_shell = 0,
1574 const unsigned int n_rotations,
1613 template <
int dim,
int spacedim>
1617 const std::string & grid_generator_function_name,
1618 const std::string & grid_generator_function_arguments);
1687 template <
int dim,
int spacedim>
1692 const double duplicated_vertex_tolerance = 1.0e-12,
1693 const bool copy_manifold_ids =
false);
1709 template <
int dim,
int spacedim>
1714 const double duplicated_vertex_tolerance = 1.0e-12,
1715 const bool copy_manifold_ids =
false);
1764 template <
int dim,
int spacedim = dim>
1767 const std::vector<unsigned int> & extents,
1806 template <
int dim,
int spacedim>
1849 template <
int dim,
int spacedim>
1916 const unsigned int n_slices,
1917 const double height,
1919 const bool copy_manifold_ids =
false,
1920 const std::vector<types::manifold_id> &manifold_priorities = {});
1931 const unsigned int n_slices,
1932 const double height,
1934 const bool copy_manifold_ids =
false,
1935 const std::vector<types::manifold_id> &manifold_priorities = {});
1959 const std::vector<double> & slice_coordinates,
1961 const bool copy_manifold_ids =
false,
1962 const std::vector<types::manifold_id> &manifold_priorities = {});
1973 const std::vector<double> & slice_coordinates,
1975 const bool copy_manifold_ids =
false,
1976 const std::vector<types::manifold_id> &manifold_priorities = {});
2012 template <
int dim,
int spacedim1,
int spacedim2>
2046 template <
int dim,
int spacedim>
2054 template <
int spacedim>
2262 template <
int dim,
int spacedim>
2266 const std::vector<unsigned int> &repetitions,
2285 template <
int dim,
int spacedim>
2288 const unsigned int repetitions,
2289 const double p1 = 0.0,
2290 const double p2 = 1.0,
2309 template <
template <
int,
int>
class MeshType,
int dim,
int spacedim>
2310 struct ExtractBoundaryMesh
2313 std::map<
typename MeshType<dim - 1, spacedim>::cell_iterator,
2314 typename MeshType<dim, spacedim>::face_iterator>;
2393 template <
template <
int,
int>
class MeshType,
int dim,
int spacedim>
2395 std::map<
typename MeshType<dim - 1, spacedim>::cell_iterator,
2396 typename MeshType<dim, spacedim>::face_iterator>
2398 typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
2401 MeshType<dim - 1, spacedim> & surface_mesh,
2402 const std::set<types::boundary_id> &boundary_ids =
2403 std::set<types::boundary_id>());
2423 <<
"The number of repetitions " << arg1 <<
" must be >=1.");
2429 <<
"The vector of repetitions must have " << arg1
2436 "The input to this function is oriented in a way that will"
2437 " cause all cells to have negative measure.");
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
__global__ void set(Number *val, const Number s, const size_type N)
#define DeclException0(Exception0)
static ::ExceptionBase & ExcInvalidInputOrientation()
static ::ExceptionBase & ExcInvalidRepetitions(int arg1)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidRadii()
#define DeclException1(Exception1, type1, outsequence)
static ::ExceptionBase & ExcInvalidRepetitionsDimension(int arg1)
void subdivided_hyper_cube_with_simplices(Triangulation< dim, spacedim > &tria, const unsigned int repetitions, const double p1=0.0, const double p2=1.0, const bool colorize=false)
void subdivided_hyper_rectangle_with_simplices(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &repetitions, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false)
void create_triangulation(Triangulation< dim, dim > &tria, const AdditionalData &additional_data=AdditionalData())
void parallelepiped(Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false)
void hyper_cross(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &sizes, const bool colorize_cells=false)
A center cell with stacks of cell protruding from each surface.
std::map< typename MeshType< dim - 1, spacedim >::cell_iterator, typename MeshType< dim, spacedim >::face_iterator > extract_boundary_mesh(const MeshType< dim, spacedim > &volume_mesh, MeshType< dim - 1, spacedim > &surface_mesh, const std::set< types::boundary_id > &boundary_ids=std::set< types::boundary_id >())
void hyper_ball_balanced(Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.)
void plate_with_a_hole(Triangulation< dim > &tria, const double inner_radius=0.4, const double outer_radius=1., const double pad_bottom=2., const double pad_top=2., const double pad_left=1., const double pad_right=1., const Point< dim > ¢er=Point< dim >(), const types::manifold_id polar_manifold_id=0, const types::manifold_id tfi_manifold_id=1, const double L=1., const unsigned int n_slices=2, const bool colorize=false)
Rectangular plate with an (offset) cylindrical hole.
void generate_from_name_and_arguments(Triangulation< dim, spacedim > &tria, const std::string &grid_generator_function_name, const std::string &grid_generator_function_arguments)
void enclosed_hyper_cube(Triangulation< dim > &tria, const double left=0., const double right=1., const double thickness=1., const bool colorize=false)
void replicate_triangulation(const Triangulation< dim, spacedim > &input, const std::vector< unsigned int > &extents, Triangulation< dim, spacedim > &result)
Replicate a given triangulation in multiple coordinate axes.
void parallelogram(Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false)
void general_cell(Triangulation< dim, spacedim > &tria, const std::vector< Point< spacedim > > &vertices, const bool colorize=false)
void subdivided_hyper_cube(Triangulation< dim, spacedim > &tria, const unsigned int repetitions, const double left=0., const double right=1., const bool colorize=false)
void hyper_L(Triangulation< dim > &tria, const double left=-1., const double right=1., const bool colorize=false)
void hyper_cube_slit(Triangulation< dim > &tria, const double left=0., const double right=1., const bool colorize=false)
void hyper_ball(Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false)
void eccentric_hyper_shell(Triangulation< dim > &triangulation, const Point< dim > &inner_center, const Point< dim > &outer_center, const double inner_radius, const double outer_radius, const unsigned int n_cells)
void hyper_rectangle(Triangulation< dim, spacedim > &tria, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false)
void cylinder(Triangulation< dim > &tria, const double radius=1., const double half_length=1.)
void moebius(Triangulation< 3, 3 > &tria, const unsigned int n_cells, const unsigned int n_rotations, const double R, const double r)
void extrude_triangulation(const Triangulation< 2, 2 > &input, const unsigned int n_slices, const double height, Triangulation< 3, 3 > &result, const bool copy_manifold_ids=false, const std::vector< types::manifold_id > &manifold_priorities={})
void half_hyper_shell(Triangulation< dim > &tria, const Point< dim > ¢er, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, const bool colorize=false)
void quarter_hyper_ball(Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.)
void cylinder_shell(Triangulation< dim > &tria, const double length, const double inner_radius, const double outer_radius, const unsigned int n_radial_cells=0, const unsigned int n_axial_cells=0)
void cheese(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &holes)
Rectangular domain with rectangular pattern of holes.
void merge_triangulations(const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result, const double duplicated_vertex_tolerance=1.0e-12, const bool copy_manifold_ids=false)
void non_standard_orientation_mesh(Triangulation< 2 > &tria, const bool rotate_left_square, const bool rotate_right_square)
void torus(Triangulation< dim, spacedim > &tria, const double R, const double r, const unsigned int n_cells_toroidal=6, const double phi=2.0 *numbers::PI)
void create_union_triangulation(const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result)
void subdivided_parallelepiped(Triangulation< dim > &tria, const unsigned int n_subdivisions, const Point< dim >(&corners)[dim], const bool colorize=false)
void subdivided_cylinder(Triangulation< dim > &tria, const unsigned int x_subdivisions, const double radius=1., const double half_length=1.)
void channel_with_cylinder(Triangulation< dim > &tria, const double shell_region_width=0.03, const unsigned int n_shells=2, const double skewness=2.0, const bool colorize=false)
void subdivided_hyper_L(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &repetitions, const Point< dim > &bottom_left, const Point< dim > &top_right, const std::vector< int > &n_cells_to_remove)
void hyper_sphere(Triangulation< spacedim - 1, spacedim > &tria, const Point< spacedim > ¢er=Point< spacedim >(), const double radius=1.)
void concentric_hyper_shells(Triangulation< dim > &triangulation, const Point< dim > ¢er, const double inner_radius=0.125, const double outer_radius=0.25, const unsigned int n_shells=1, const double skewness=0.1, const unsigned int n_cells_per_shell=0, const bool colorize=false)
void convert_hypercube_to_simplex_mesh(const Triangulation< dim, spacedim > &in_tria, Triangulation< dim, spacedim > &out_tria)
void subdivided_hyper_rectangle(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &repetitions, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false)
void quarter_hyper_shell(Triangulation< dim > &tria, const Point< dim > ¢er, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, const bool colorize=false)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
void hyper_shell(Triangulation< dim > &tria, const Point< dim > ¢er, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, bool colorize=false)
void create_triangulation_with_removed_cells(const Triangulation< dim, spacedim > &input_triangulation, const std::set< typename Triangulation< dim, spacedim >::active_cell_iterator > &cells_to_remove, Triangulation< dim, spacedim > &result)
void simplex(Triangulation< dim, dim > &tria, const std::vector< Point< dim > > &vertices)
void hyper_cube_with_cylindrical_hole(Triangulation< dim > &triangulation, const double inner_radius=.25, const double outer_radius=.5, const double L=.5, const unsigned int repetitions=1, const bool colorize=false)
void truncated_cone(Triangulation< dim > &tria, const double radius_0=1.0, const double radius_1=0.5, const double half_length=1.0)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
void half_hyper_ball(Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.)
void flatten_triangulation(const Triangulation< dim, spacedim1 > &in_tria, Triangulation< dim, spacedim2 > &out_tria)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
static constexpr double PI
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
unsigned int n_subdivision_x_0
unsigned int n_subdivision_y
unsigned int n_subdivision_x_2
void add_parameters(ParameterHandler &prm)
unsigned int airfoil_sampling_factor
unsigned int n_subdivision_x_1
Point< 2, double > joukowski_center