16#ifndef dealii_grid_generator_h
17#define dealii_grid_generator_h
95 template <
int dim,
int spacedim>
98 const double left = 0.,
99 const double right = 1.,
139 template <
int dim,
int spacedim>
171 template <
int dim,
int spacedim>
174 const unsigned int repetitions,
175 const double left = 0.,
176 const double right = 1.,
214 template <
int dim,
int spacedim>
270 template <
int dim,
int spacedim>
273 const std::vector<unsigned int> &repetitions,
296 const std::vector<std::vector<double>> &step_sizes,
314 const std::vector<std::vector<double>> &spacing,
343 template <
int dim,
int spacedim>
346 const std::vector<unsigned int> &holes);
394 const double inner_radius = 0.4,
395 const double outer_radius = 1.,
396 const double pad_bottom = 2.,
397 const double pad_top = 2.,
398 const double pad_left = 1.,
399 const double pad_right = 1.,
404 const unsigned int n_slices = 2,
496 const double shell_region_width = 0.03,
497 const unsigned int n_shells = 2,
498 const double skewness = 2.0,
521 template <
int dim,
int spacedim>
584 const unsigned int n_subdivisions,
600 const unsigned int (&n_subdivisions)[dim],
602 const unsigned int *n_subdivisions,
630 template <
int dim,
int spacedim>
635 const std::vector<unsigned int> &subdivisions = {},
657 const double left = 0.,
658 const double right = 1.,
659 const double thickness = 1.,
742 const double radius = 1.,
743 const bool attach_spherical_manifold_on_boundary_cells =
false);
781 const double radius = 1.);
803 const unsigned int n_rotate_middle_square);
828 const bool face_orientation,
829 const bool face_flip,
830 const bool face_rotation,
831 const bool manipulate_left_cube);
859 template <
int spacedim>
863 const double radius = 1.);
902 const double radius = 1.);
920 const double radius = 1.);
948 const double radius = 1.,
949 const double half_length = 1.);
988 const unsigned int x_subdivisions,
989 const double radius = 1.,
990 const double half_length = 1.);
1026 const double radius_0 = 1.0,
1027 const double radius_1 = 0.5,
1028 const double half_length = 1.0);
1153 template <
int dim,
int spacedim>
1158 const double aspect_ratio = 0.5);
1186 template <
int dim,
int spacedim>
1189 const std::vector<unsigned int> &sizes,
1190 const bool colorize_cells =
false);
1230 const double left = -1.,
1231 const double right = 1.,
1265 template <
int dim,
int spacedim>
1268 const std::vector<unsigned int> &repetitions,
1271 const std::vector<int> & n_cells_to_remove);
1295 const double left = 0.,
1296 const double right = 1.,
1360 const double inner_radius,
1361 const double outer_radius,
1362 const unsigned int n_cells = 0,
1403 const double inner_radius,
1404 const double outer_radius,
1405 const unsigned int n_cells);
1440 const double inner_radius,
1441 const double outer_radius,
1442 const unsigned int n_cells = 0,
1475 const double inner_radius,
1476 const double outer_radius,
1477 const unsigned int n_cells = 0,
1509 const double length,
1510 const double inner_radius,
1511 const double outer_radius,
1512 const unsigned int n_radial_cells = 0,
1513 const unsigned int n_axial_cells = 0);
1564 template <
int dim,
int spacedim>
1569 const unsigned int n_cells_toroidal = 6,
1604 const double inner_radius = .25,
1605 const double outer_radius = .5,
1606 const double L = .5,
1607 const unsigned int repetitions = 1,
1684 const double inner_radius = 0.125,
1685 const double outer_radius = 0.25,
1686 const unsigned int n_shells = 1,
1687 const double skewness = 0.1,
1688 const unsigned int n_cells_per_shell = 0,
1706 const unsigned int n_cells,
1707 const unsigned int n_rotations,
1746 template <
int dim,
int spacedim>
1750 const std::string & grid_generator_function_name,
1751 const std::string & grid_generator_function_arguments);
1809 const double & outer_ball_radius = 1.0);
1900 template <
int dim,
int spacedim>
1905 const double duplicated_vertex_tolerance = 1.0e-12,
1906 const bool copy_manifold_ids =
false);
1922 template <
int dim,
int spacedim>
1927 const double duplicated_vertex_tolerance = 1.0e-12,
1928 const bool copy_manifold_ids =
false);
1977 template <
int dim,
int spacedim = dim>
1980 const std::vector<unsigned int> & extents,
2020 template <
int dim,
int spacedim>
2063 template <
int dim,
int spacedim>
2130 const unsigned int n_slices,
2131 const double height,
2133 const bool copy_manifold_ids =
false,
2134 const std::vector<types::manifold_id> &manifold_priorities = {});
2145 const unsigned int n_slices,
2146 const double height,
2148 const bool copy_manifold_ids =
false,
2149 const std::vector<types::manifold_id> &manifold_priorities = {});
2173 const std::vector<double> & slice_coordinates,
2175 const bool copy_manifold_ids =
false,
2176 const std::vector<types::manifold_id> &manifold_priorities = {});
2187 const std::vector<double> & slice_coordinates,
2189 const bool copy_manifold_ids =
false,
2190 const std::vector<types::manifold_id> &manifold_priorities = {});
2230 template <
int dim,
int spacedim1,
int spacedim2>
2268 template <
int dim,
int spacedim>
2276 template <
int spacedim>
2484 template <
int dim,
int spacedim>
2486 subdivided_hyper_rectangle_with_simplices(
2488 const std::vector<unsigned int> &repetitions,
2507 template <
int dim,
int spacedim>
2510 const unsigned int repetitions,
2511 const double p1 = 0.0,
2512 const double p2 = 1.0,
2531 template <
template <
int,
int>
class MeshType,
int dim,
int spacedim>
2532 struct ExtractBoundaryMesh
2535 std::map<
typename MeshType<dim - 1, spacedim>::cell_iterator,
2536 typename MeshType<dim, spacedim>::face_iterator>;
2634 template <
template <
int,
int>
class MeshType,
int dim,
int spacedim>
2639 std::map<
typename MeshType<dim - 1, spacedim>::cell_iterator,
2640 typename MeshType<dim, spacedim>::face_iterator>
2642 typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
2646 MeshType<dim - 1, spacedim> & surface_mesh,
2647 const std::set<types::boundary_id> &boundary_ids =
2648 std::set<types::boundary_id>());
2668 <<
"The number of repetitions " << arg1 <<
" must be >=1.");
2674 <<
"The vector of repetitions must have " << arg1
2681 "The input to this function is oriented in a way that will"
2682 " cause all cells to have negative measure.");
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#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 implicit_function(Triangulation< dim, 3 > &tria, const Function< 3 > &implicit_function, const CGALWrappers::AdditionalData< dim > &data=CGALWrappers::AdditionalData< dim >{}, const Point< 3 > &interior_point=Point< 3 >(), const double &outer_ball_radius=1.0)
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 convert_hypercube_to_simplex_mesh(const Triangulation< 1, spacedim > &in_tria, Triangulation< 1, spacedim > &out_tria)
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.
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 surface_mesh_to_volumetric_mesh(const Triangulation< 2, 3 > &surface_tria, Triangulation< 3 > &vol_tria, const CGALWrappers::AdditionalData< 3 > &data=CGALWrappers::AdditionalData< 3 >{})
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 pipe_junction(Triangulation< dim, spacedim > &tria, const std::vector< std::pair< Point< spacedim >, double > > &openings, const std::pair< Point< spacedim >, double > &bifurcation, const double aspect_ratio=0.5)
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 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 non_standard_orientation_mesh(Triangulation< 2 > &tria, const unsigned int n_rotate_middle_square)
return_type 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 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 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)
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
const ::Triangulation< dim, spacedim > & tria