Reference documentation for deal.II version 9.3.0
\(\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\}}\)
grid_generator.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 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_generator_h
17 #define dealii_grid_generator_h
18 
19 
20 #include <deal.II/base/config.h>
21 
23 #include <deal.II/base/function.h>
25 #include <deal.II/base/point.h>
26 #include <deal.II/base/table.h>
27 
29 #include <deal.II/grid/tria.h>
30 
31 #include <array>
32 #include <map>
33 
35 
58 namespace GridGenerator
59 {
63 
93  template <int dim, int spacedim>
94  void
96  const double left = 0.,
97  const double right = 1.,
98  const bool colorize = false);
99 
127  template <int dim>
128  void
130  const std::vector<Point<dim>> &vertices);
131 
132  /*
133  * Create a (coarse) grid with a single cell of the shape of the provided
134  * reference cell. This is a generalization of the hyper_cube() and simplex()
135  * functions above.
136  */
137  template <int dim, int spacedim>
138  void
140  const ReferenceCell & reference_cell);
141 
142 
169  template <int dim, int spacedim>
170  void
172  const unsigned int repetitions,
173  const double left = 0.,
174  const double right = 1.,
175  const bool colorize = false);
176 
212  template <int dim, int spacedim>
213  void
215  const Point<dim> & p1,
216  const Point<dim> & p2,
217  const bool colorize = false);
218 
268  template <int dim, int spacedim>
269  void
271  const std::vector<unsigned int> &repetitions,
272  const Point<dim> & p1,
273  const Point<dim> & p2,
274  const bool colorize = false);
275 
291  template <int dim>
292  void
294  const std::vector<std::vector<double>> &step_sizes,
295  const Point<dim> & p_1,
296  const Point<dim> & p_2,
297  const bool colorize = false);
298 
309  template <int dim>
310  void
312  const std::vector<std::vector<double>> &spacing,
313  const Point<dim> & p,
315  const bool colorize = false);
316 
341  template <int dim, int spacedim>
342  void
344  const std::vector<unsigned int> &holes);
345 
389  template <int dim>
390  void
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.,
398  const Point<dim> & center = Point<dim>(),
399  const types::manifold_id polar_manifold_id = 0,
400  const types::manifold_id tfi_manifold_id = 1,
401  const double L = 1.,
402  const unsigned int n_slices = 2,
403  const bool colorize = false);
404 
491  template <int dim>
492  void
494  const double shell_region_width = 0.03,
495  const unsigned int n_shells = 2,
496  const double skewness = 2.0,
497  const bool colorize = false);
498 
519  template <int dim, int spacedim>
520  void
522  const std::vector<Point<spacedim>> &vertices,
523  const bool colorize = false);
524 
541  template <int dim>
542  void
544  const Point<dim> (&corners)[dim],
545  const bool colorize = false);
546 
562  template <int dim>
563  void
565  const Point<dim> (&corners)[dim],
566  const bool colorize = false);
567 
579  template <int dim>
580  void
582  const unsigned int n_subdivisions,
583  const Point<dim> (&corners)[dim],
584  const bool colorize = false);
585 
594  template <int dim>
595  void
597 #ifndef _MSC_VER
598  const unsigned int (&n_subdivisions)[dim],
599 #else
600  const unsigned int *n_subdivisions,
601 #endif
602  const Point<dim> (&corners)[dim],
603  const bool colorize = false);
604 
628  template <int dim, int spacedim>
629  void
631  const Point<spacedim> & origin,
632  const std::array<Tensor<1, spacedim>, dim> &edges,
633  const std::vector<unsigned int> &subdivisions = {},
634  const bool colorize = false);
635 
652  template <int dim>
653  void
655  const double left = 0.,
656  const double right = 1.,
657  const double thickness = 1.,
658  const bool colorize = false);
659 
736  template <int dim>
737  void
739  const Point<dim> & center = Point<dim>(),
740  const double radius = 1.,
741  const bool attach_spherical_manifold_on_boundary_cells = false);
742 
775  template <int dim>
776  void
778  const Point<dim> & center = Point<dim>(),
779  const double radius = 1.);
780 
802  const bool rotate_left_square,
803  const bool rotate_right_square);
804 
827  const bool face_orientation,
828  const bool face_flip,
829  const bool face_rotation,
830  const bool manipulate_left_cube);
831 
832 
858  template <int spacedim>
861  const double radius = 1.);
862 
896  template <int dim>
897  void
899  const Point<dim> & center = Point<dim>(),
900  const double radius = 1.);
901 
914  template <int dim>
915  void
917  const Point<dim> & center = Point<dim>(),
918  const double radius = 1.);
919 
943  template <int dim>
944  void
946  const double radius = 1.,
947  const double half_length = 1.);
948 
949 
983  template <int dim>
984  void
986  const unsigned int x_subdivisions,
987  const double radius = 1.,
988  const double half_length = 1.);
989 
990 
1021  template <int dim>
1022  void
1024  const double radius_0 = 1.0,
1025  const double radius_1 = 0.5,
1026  const double half_length = 1.0);
1027 
1054  template <int dim, int spacedim>
1055  void
1057  const std::vector<unsigned int> &sizes,
1058  const bool colorize_cells = false);
1059 
1095  template <int dim>
1096  void
1098  const double left = -1.,
1099  const double right = 1.,
1100  const bool colorize = false);
1101 
1133  template <int dim, int spacedim>
1134  void
1136  const std::vector<unsigned int> &repetitions,
1137  const Point<dim> & bottom_left,
1138  const Point<dim> & top_right,
1139  const std::vector<int> & n_cells_to_remove);
1140 
1160  template <int dim>
1161  void
1163  const double left = 0.,
1164  const double right = 1.,
1165  const bool colorize = false);
1166 
1224  template <int dim>
1225  void
1227  const Point<dim> & center,
1228  const double inner_radius,
1229  const double outer_radius,
1230  const unsigned int n_cells = 0,
1231  bool colorize = false);
1232 
1266  template <int dim>
1267  void
1269  const Point<dim> & inner_center,
1270  const Point<dim> & outer_center,
1271  const double inner_radius,
1272  const double outer_radius,
1273  const unsigned int n_cells);
1274 
1304  template <int dim>
1305  void
1307  const Point<dim> & center,
1308  const double inner_radius,
1309  const double outer_radius,
1310  const unsigned int n_cells = 0,
1311  const bool colorize = false);
1312 
1313 
1339  template <int dim>
1340  void
1342  const Point<dim> & center,
1343  const double inner_radius,
1344  const double outer_radius,
1345  const unsigned int n_cells = 0,
1346  const bool colorize = false);
1347 
1374  template <int dim>
1375  void
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);
1382 
1432  template <int dim, int spacedim>
1433  void
1435  const double R,
1436  const double r,
1437  const unsigned int n_cells_toroidal = 6,
1438  const double phi = 2.0 * numbers::PI);
1439 
1469  template <int dim>
1470  void
1472  const double inner_radius = .25,
1473  const double outer_radius = .5,
1474  const double L = .5,
1475  const unsigned int repetitions = 1,
1476  const bool colorize = false);
1477 
1548  template <int dim>
1549  void
1551  const Point<dim> & center,
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,
1557  const bool colorize = false);
1558 
1572  void moebius(Triangulation<3, 3> &tria,
1573  const unsigned int n_cells,
1574  const unsigned int n_rotations,
1575  const double R,
1576  const double r);
1577 
1613  template <int dim, int spacedim>
1614  void
1617  const std::string & grid_generator_function_name,
1618  const std::string & grid_generator_function_arguments);
1620 
1624 
1687  template <int dim, int spacedim>
1688  void
1689  merge_triangulations(const Triangulation<dim, spacedim> &triangulation_1,
1690  const Triangulation<dim, spacedim> &triangulation_2,
1692  const double duplicated_vertex_tolerance = 1.0e-12,
1693  const bool copy_manifold_ids = false);
1694 
1709  template <int dim, int spacedim>
1710  void
1712  const std::vector<const Triangulation<dim, spacedim> *> &triangulations,
1714  const double duplicated_vertex_tolerance = 1.0e-12,
1715  const bool copy_manifold_ids = false);
1716 
1764  template <int dim, int spacedim = dim>
1765  void
1767  const std::vector<unsigned int> & extents,
1768  Triangulation<dim, spacedim> & result);
1769 
1806  template <int dim, int spacedim>
1807  void
1809  const Triangulation<dim, spacedim> &triangulation_1,
1810  const Triangulation<dim, spacedim> &triangulation_2,
1811  Triangulation<dim, spacedim> & result);
1812 
1849  template <int dim, int spacedim>
1850  void
1852  const Triangulation<dim, spacedim> &input_triangulation,
1854  & cells_to_remove,
1856 
1913  void
1915  const Triangulation<2, 2> & input,
1916  const unsigned int n_slices,
1917  const double height,
1918  Triangulation<3, 3> & result,
1919  const bool copy_manifold_ids = false,
1920  const std::vector<types::manifold_id> &manifold_priorities = {});
1921 
1928  void
1930  const Triangulation<2, 2> & input,
1931  const unsigned int n_slices,
1932  const double height,
1933  Triangulation<2, 2> & result,
1934  const bool copy_manifold_ids = false,
1935  const std::vector<types::manifold_id> &manifold_priorities = {});
1936 
1937 
1956  void
1958  const Triangulation<2, 2> & input,
1959  const std::vector<double> & slice_coordinates,
1960  Triangulation<3, 3> & result,
1961  const bool copy_manifold_ids = false,
1962  const std::vector<types::manifold_id> &manifold_priorities = {});
1963 
1970  void
1972  const Triangulation<2, 2> & input,
1973  const std::vector<double> & slice_coordinates,
1974  Triangulation<2, 2> & result,
1975  const bool copy_manifold_ids = false,
1976  const std::vector<types::manifold_id> &manifold_priorities = {});
1977 
1978 
1979 
2012  template <int dim, int spacedim1, int spacedim2>
2013  void
2015  Triangulation<dim, spacedim2> & out_tria);
2016 
2046  template <int dim, int spacedim>
2047  void
2049  Triangulation<dim, spacedim> &out_tria);
2050 
2054  template <int spacedim>
2055  void
2057  Triangulation<1, spacedim> & out_tria);
2058 
2059 
2066  namespace Airfoil
2067  {
2073  {
2078  std::string airfoil_type;
2079 
2087  std::string naca_id;
2088 
2095 
2101 
2106  double height;
2107 
2111  double length_b2;
2112 
2129 
2136  double bias_factor;
2137 
2141  unsigned int refinements;
2142 
2146  unsigned int n_subdivision_x_0;
2147 
2151  unsigned int n_subdivision_x_1;
2152 
2156  unsigned int n_subdivision_x_2;
2157 
2161  unsigned int n_subdivision_y;
2162 
2175 
2179  AdditionalData();
2180 
2186  void
2188  };
2189 
2219  template <int dim>
2220  void
2223  const AdditionalData & additional_data = AdditionalData());
2224 
2225 
2226 
2239  template <int dim>
2240  void
2242  Triangulation<dim, dim> & tria,
2243  std::vector<GridTools::PeriodicFacePair<
2244  typename Triangulation<dim, dim>::cell_iterator>> &periodic_faces,
2245  const AdditionalData &additional_data = AdditionalData());
2246 
2247  } // namespace Airfoil
2248 
2262  template <int dim, int spacedim>
2263  void
2266  const std::vector<unsigned int> &repetitions,
2267  const Point<dim> & p1,
2268  const Point<dim> & p2,
2269  const bool colorize = false);
2270 
2285  template <int dim, int spacedim>
2286  void
2288  const unsigned int repetitions,
2289  const double p1 = 0.0,
2290  const double p2 = 1.0,
2291  const bool colorize = false);
2292 
2294 
2300 
2302 #ifdef _MSC_VER
2303  // Microsoft's VC++ has a bug where it doesn't want to recognize that
2304  // an implementation (definition) of the extract_boundary_mesh function
2305  // matches a declaration. This can apparently only be avoided by
2306  // doing some contortion with the return type using the following
2307  // intermediate type. This is only used when using MS VC++ and uses
2308  // the direct way of doing it otherwise
2309  template <template <int, int> class MeshType, int dim, int spacedim>
2310  struct ExtractBoundaryMesh
2311  {
2312  using return_type =
2313  std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2314  typename MeshType<dim, spacedim>::face_iterator>;
2315  };
2316 #endif
2317 
2393  template <template <int, int> class MeshType, int dim, int spacedim>
2394 #ifndef _MSC_VER
2395  std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2396  typename MeshType<dim, spacedim>::face_iterator>
2397 #else
2398  typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
2399 #endif
2400  extract_boundary_mesh(const MeshType<dim, spacedim> & volume_mesh,
2401  MeshType<dim - 1, spacedim> & surface_mesh,
2402  const std::set<types::boundary_id> &boundary_ids =
2403  std::set<types::boundary_id>());
2404 
2406 
2407 
2411 
2413 
2422  int,
2423  << "The number of repetitions " << arg1 << " must be >=1.");
2428  int,
2429  << "The vector of repetitions must have " << arg1
2430  << " elements.");
2431 
2436  "The input to this function is oriented in a way that will"
2437  " cause all cells to have negative measure.");
2439 
2440 #ifndef DOXYGEN
2441  // These functions are only implemented with specializations; declare them
2442  // here
2443  template <>
2445  const double,
2446  const double,
2447  const double,
2448  const unsigned int,
2449  const bool);
2450 
2451  template <>
2453  const double,
2454  const double,
2455  const double,
2456  const unsigned int,
2457  const bool);
2458 
2459  template <>
2461  const double,
2462  const double,
2463  const double,
2464  const unsigned int,
2465  const bool);
2466 
2467  template <>
2469  const double,
2470  const unsigned int,
2471  const double,
2472  const bool);
2473 
2474  template <>
2476  const double,
2477  const unsigned int,
2478  const double,
2479  const bool);
2480 
2481  template <>
2483  const double,
2484  const unsigned int,
2485  const double,
2486  const bool);
2487 #endif
2488 } // namespace GridGenerator
2489 
2490 
2491 
2493 
2494 #endif
void create_union_triangulation(const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result)
static const char L
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 simplex(Triangulation< dim, dim > &tria, const std::vector< Point< dim >> &vertices)
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)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
void concentric_hyper_shells(Triangulation< dim > &triangulation, const Point< dim > &center, 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 general_cell(Triangulation< dim, spacedim > &tria, const std::vector< Point< spacedim >> &vertices, 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 hyper_rectangle(Triangulation< dim, spacedim > &tria, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false)
void moebius(Triangulation< 3, 3 > &tria, const unsigned int n_cells, const unsigned int n_rotations, const double R, const double r)
void flatten_triangulation(const Triangulation< dim, spacedim1 > &in_tria, Triangulation< dim, spacedim2 > &out_tria)
void quarter_hyper_shell(Triangulation< dim > &tria, const Point< dim > &center, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, const bool colorize=false)
void enclosed_hyper_cube(Triangulation< dim > &tria, const double left=0., const double right=1., const double thickness=1., const bool colorize=false)
void subdivided_parallelepiped(Triangulation< dim > &tria, const unsigned int n_subdivisions, const Point< dim >(&corners)[dim], const bool colorize=false)
static ::ExceptionBase & ExcInvalidRepetitionsDimension(int arg1)
unsigned int material_id
Definition: types.h:152
void parallelogram(Triangulation< dim > &tria, const Point< dim >(&corners)[dim], 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 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 > &center=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 hyper_shell(Triangulation< dim > &tria, const Point< dim > &center, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, bool colorize=false)
static ::ExceptionBase & ExcInvalidInputOrientation()
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 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 create_triangulation(Triangulation< dim, dim > &tria, const AdditionalData &additional_data=AdditionalData())
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_L(Triangulation< dim > &tria, const double left=-1., const double right=1., const bool colorize=false)
void half_hyper_ball(Triangulation< dim > &tria, const Point< dim > &center=Point< dim >(), const double radius=1.)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
void hyper_sphere(Triangulation< spacedim - 1, spacedim > &tria, const Point< spacedim > &center=Point< spacedim >(), const double radius=1.)
void non_standard_orientation_mesh(Triangulation< 2 > &tria, const bool rotate_left_square, const bool rotate_right_square)
void cylinder(Triangulation< dim > &tria, const double radius=1., const double half_length=1.)
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
void hyper_ball_balanced(Triangulation< dim > &tria, const Point< dim > &center=Point< dim >(), const double radius=1.)
#define DeclException0(Exception0)
Definition: exceptions.h:470
void quarter_hyper_ball(Triangulation< dim > &tria, const Point< dim > &center=Point< dim >(), const double radius=1.)
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:395
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 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 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={})
Point< 3 > vertices[4]
void add_parameters(ParameterHandler &prm)
void hyper_cube_slit(Triangulation< dim > &tria, const double left=0., const double right=1., const bool colorize=false)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
Definition: tria.cc:12587
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 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 reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
void generate_from_name_and_arguments(Triangulation< dim, spacedim > &tria, const std::string &grid_generator_function_name, const std::string &grid_generator_function_arguments)
static ::ExceptionBase & ExcInvalidRadii()
void cheese(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &holes)
Rectangular domain with rectangular pattern of holes.
__global__ void set(Number *val, const Number s, const size_type N)
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_cylinder(Triangulation< dim > &tria, const unsigned int x_subdivisions, const double radius=1., const double half_length=1.)
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 parallelepiped(Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false)
void half_hyper_shell(Triangulation< dim > &tria, const Point< dim > &center, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, const bool colorize=false)
static constexpr double PI
Definition: numbers.h:231
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:394
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 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 hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
Point< 3 > center
void hyper_ball(Triangulation< dim > &tria, const Point< dim > &center=Point< dim >(), const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false)
bool colorize
Definition: grid_out.cc:4589
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
void convert_hypercube_to_simplex_mesh(const Triangulation< dim, spacedim > &in_tria, Triangulation< dim, spacedim > &out_tria)
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)
static ::ExceptionBase & ExcInvalidRepetitions(int arg1)