Reference documentation for deal.II version 9.1.0-pre
grid_generator.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1999 - 2018 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 
22 #include <deal.II/base/exceptions.h>
23 #include <deal.II/base/function.h>
24 #include <deal.II/base/point.h>
25 #include <deal.II/base/table.h>
26 
27 #include <deal.II/grid/tria.h>
28 
29 #include <array>
30 #include <map>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
56 namespace GridGenerator
57 {
61 
89  template <int dim, int spacedim>
90  void
92  const double left = 0.,
93  const double right = 1.,
94  const bool colorize = false);
95 
120  template <int dim>
121  void
123  const std::vector<Point<dim>> &vertices);
124 
138  template <int dim, int spacedim>
139  void
141  const unsigned int repetitions,
142  const double left = 0.,
143  const double right = 1.);
144 
176  template <int dim, int spacedim>
177  void
179  const Point<dim> & p1,
180  const Point<dim> & p2,
181  const bool colorize = false);
182 
232  template <int dim, int spacedim>
233  void
235  const std::vector<unsigned int> &repetitions,
236  const Point<dim> & p1,
237  const Point<dim> & p2,
238  const bool colorize = false);
239 
255  template <int dim>
256  void
258  const std::vector<std::vector<double>> &step_sizes,
259  const Point<dim> & p_1,
260  const Point<dim> & p_2,
261  const bool colorize = false);
262 
273  template <int dim>
274  void
276  const std::vector<std::vector<double>> &spacing,
277  const Point<dim> & p,
278  const Table<dim, types::material_id> &material_id,
279  const bool colorize = false);
280 
307  template <int dim, int spacedim>
308  void
310  const std::vector<unsigned int> &holes);
311 
354  template <int dim>
355  void
357  const double inner_radius = 0.4,
358  const double outer_radius = 1.,
359  const double pad_bottom = 2.,
360  const double pad_top = 2.,
361  const double pad_left = 1.,
362  const double pad_right = 1.,
363  const Point<dim> center = Point<dim>(),
364  const types::manifold_id polar_manifold_id = 0,
365  const types::manifold_id tfi_manifold_id = 1,
366  const double L = 1.,
367  const unsigned int n_slices = 2,
368  const bool colorize = false);
369 
385  template <int dim>
386  void
388  const std::vector<Point<dim>> &vertices,
389  const bool colorize = false);
390 
401  template <int dim>
402  void
404  const Point<dim> (&corners)[dim],
405  const bool colorize = false);
406 
421  template <int dim>
422  void
424  const Point<dim> (&corners)[dim],
425  const bool colorize = false);
426 
437  template <int dim>
438  void
440  const unsigned int n_subdivisions,
441  const Point<dim> (&corners)[dim],
442  const bool colorize = false);
443 
451  template <int dim>
452  void
454 #ifndef _MSC_VER
455  const unsigned int (&n_subdivisions)[dim],
456 #else
457  const unsigned int *n_subdivisions,
458 #endif
459  const Point<dim> (&corners)[dim],
460  const bool colorize = false);
461 
485  template <int dim, int spacedim>
486  void
489  const Point<spacedim> & origin,
490  const std::array<Tensor<1, spacedim>, dim> &edges,
491  const std::vector<unsigned int> &subdivisions = std::vector<unsigned int>(),
492  const bool colorize = false);
493 
509  template <int dim>
510  void
512  const double left = 0.,
513  const double right = 1.,
514  const double thickness = 1.,
515  const bool colorize = false);
516 
551  template <int dim>
552  void
554  const Point<dim> & center = Point<dim>(),
555  const double radius = 1.,
556  const bool attach_spherical_manifold_on_boundary_cells = false);
557 
583  template <int spacedim>
585  const Point<spacedim> &center = Point<spacedim>(),
586  const double radius = 1.);
587 
599  template <int dim>
600  void
602  const Point<dim> & center = Point<dim>(),
603  const double radius = 1.);
604 
616  template <int dim>
617  void
619  const Point<dim> & center = Point<dim>(),
620  const double radius = 1.);
621 
646  template <int dim>
647  void
649  const double radius = 1.,
650  const double half_length = 1.);
651 
678  template <int dim>
679  void
681  const double radius_0 = 1.0,
682  const double radius_1 = 0.5,
683  const double half_length = 1.0);
684 
714  template <int dim, int spacedim>
715  void
717  const std::vector<unsigned int> &sizes,
718  const bool colorize_cells = false);
719 
754  template <int dim>
755  void
757  const double left = -1.,
758  const double right = 1.,
759  const bool colorize = false);
760 
779  template <int dim>
780  void
782  const double left = 0.,
783  const double right = 1.,
784  const bool colorize = false);
785 
824  template <int dim>
825  void
827  const Point<dim> & center,
828  const double inner_radius,
829  const double outer_radius,
830  const unsigned int n_cells = 0,
831  bool colorize = false);
832 
858  template <int dim>
859  void
861  const Point<dim> & center,
862  const double inner_radius,
863  const double outer_radius,
864  const unsigned int n_cells = 0,
865  const bool colorize = false);
866 
867 
892  template <int dim>
893  void
895  const Point<dim> & center,
896  const double inner_radius,
897  const double outer_radius,
898  const unsigned int n_cells = 0,
899  const bool colorize = false);
900 
920  template <int dim>
921  void
923  const double length,
924  const double inner_radius,
925  const double outer_radius,
926  const unsigned int n_radial_cells = 0,
927  const unsigned int n_axial_cells = 0);
928 
929 
930 
952  template <int dim, int spacedim>
953  void
954  torus(Triangulation<dim, spacedim> &tria, const double R, const double r);
955 
985  template <int dim>
986  void
988  const double inner_radius = .25,
989  const double outer_radius = .5,
990  const double L = .5,
991  const unsigned int repetitions = 1,
992  const bool colorize = false);
993 
1065  template <int dim>
1066  void
1068  const Point<dim> & center,
1069  const double inner_radius = 0.125,
1070  const double outer_radius = 0.25,
1071  const unsigned int n_shells = 1,
1072  const double skewness = 0.1,
1073  const unsigned int n_cells_per_shell = 0,
1074  const bool colorize = false);
1075 
1089  void moebius(Triangulation<3, 3> &tria,
1090  const unsigned int n_cells,
1091  const unsigned int n_rotations,
1092  const double R,
1093  const double r);
1094 
1096 
1100 
1159  template <int dim, int spacedim>
1160  void
1161  merge_triangulations(const Triangulation<dim, spacedim> &triangulation_1,
1162  const Triangulation<dim, spacedim> &triangulation_2,
1164  const double duplicated_vertex_tolerance = 1.0e-12,
1165  const bool copy_manifold_ids = false);
1166 
1181  template <int dim, int spacedim>
1182  void
1184  const std::initializer_list<const Triangulation<dim, spacedim> *const>
1185  triangulations,
1187  const double duplicated_vertex_tolerance = 1.0e-12,
1188  const bool copy_manifold_ids = false);
1189 
1223  template <int dim, int spacedim>
1224  void
1226  const Triangulation<dim, spacedim> &triangulation_1,
1227  const Triangulation<dim, spacedim> &triangulation_2,
1228  Triangulation<dim, spacedim> & result);
1229 
1265  template <int dim, int spacedim>
1266  void
1268  const Triangulation<dim, spacedim> &input_triangulation,
1270  & cells_to_remove,
1272 
1287  void
1289  const unsigned int n_slices,
1290  const double height,
1291  Triangulation<3, 3> & result);
1292 
1311  void
1313  const std::vector<double> &slice_coordinates,
1314  Triangulation<3, 3> & result);
1315 
1316 
1347  template <int dim, int spacedim1, int spacedim2>
1348  void
1350  Triangulation<dim, spacedim2> & out_tria);
1351 
1353 
1358 
1360 #ifdef _MSC_VER
1361  // Microsoft's VC++ has a bug where it doesn't want to recognize that
1362  // an implementation (definition) of the extract_boundary_mesh function
1363  // matches a declaration. This can apparently only be avoided by
1364  // doing some contortion with the return type using the following
1365  // intermediate type. This is only used when using MS VC++ and uses
1366  // the direct way of doing it otherwise
1367  template <template <int, int> class MeshType, int dim, int spacedim>
1368  struct ExtractBoundaryMesh
1369  {
1370  using return_type =
1371  std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
1372  typename MeshType<dim, spacedim>::face_iterator>;
1373  };
1374 #endif
1375 
1450  template <template <int, int> class MeshType, int dim, int spacedim>
1451 #ifndef _MSC_VER
1452  std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
1453  typename MeshType<dim, spacedim>::face_iterator>
1454 #else
1455  typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
1456 #endif
1457  extract_boundary_mesh(const MeshType<dim, spacedim> & volume_mesh,
1458  MeshType<dim - 1, spacedim> & surface_mesh,
1459  const std::set<types::boundary_id> &boundary_ids =
1460  std::set<types::boundary_id>());
1461 
1463 
1464 
1468 
1470 
1479  int,
1480  << "The number of repetitions " << arg1 << " must be >=1.");
1485  int,
1486  << "The vector of repetitions must have " << arg1
1487  << " elements.");
1488 
1493  "The input to this function is oriented in a way that will"
1494  " cause all cells to have negative measure.");
1496 
1497 #ifndef DOXYGEN
1498  // These functions are only implemented with specializations; declare them
1499  // here
1500  template <>
1502  const double,
1503  const double,
1504  const double,
1505  const unsigned int,
1506  const bool);
1507 
1508  template <>
1510  const double,
1511  const double,
1512  const double,
1513  const unsigned int,
1514  const bool);
1515 
1516  template <>
1518  const double,
1519  const double,
1520  const double,
1521  const unsigned int,
1522  const bool);
1523 #endif
1524 
1525 } // namespace GridGenerator
1526 
1527 
1528 
1529 DEAL_II_NAMESPACE_CLOSE
1530 
1531 #endif
void create_union_triangulation(const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result)
unsigned int manifold_id
Definition: types.h:123
void torus(Triangulation< dim, spacedim > &tria, const double R, const double r)
void simplex(Triangulation< dim, dim > &tria, const std::vector< Point< dim >> &vertices)
Triangulation of a d-simplex with (d+1) vertices and mesh cells.
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 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 extrude_triangulation(const Triangulation< 2, 2 > &input, const unsigned int n_slices, const double height, Triangulation< 3, 3 > &result)
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)
void parallelogram(Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false)
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)
void hyper_sphere(Triangulation< spacedim-1, spacedim > &tria, const Point< spacedim > &center=Point< spacedim >(), const double radius=1.)
static::ExceptionBase & ExcInvalidInputOrientation()
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 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:408
void cylinder(Triangulation< dim > &tria, const double radius=1., const double half_length=1.)
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:397
#define DeclException0(Exception0)
Definition: exceptions.h:385
void quarter_hyper_ball(Triangulation< dim > &tria, const Point< dim > &center=Point< dim >(), const double radius=1.)
void subdivided_hyper_cube(Triangulation< dim, spacedim > &tria, const unsigned int repetitions, const double left=0., const double right=1.)
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)
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_cube_slit(Triangulation< dim > &tria, const double left=0., const double right=1., 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 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)
static::ExceptionBase & ExcInvalidRadii()
void cheese(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &holes)
Rectangular domain with rectangular pattern of holes.
void parallelepiped(Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false)
void general_cell(Triangulation< dim > &tria, const std::vector< Point< dim >> &vertices, 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)
void hyper_cube(Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false)
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)
Definition: table.h:37
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_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)