Reference documentation for deal.II version 9.6.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\}}\)
Loading...
Searching...
No Matches
grid_generator.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1999 - 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_generator_h
16#define dealii_grid_generator_h
17
18
19#include <deal.II/base/config.h>
20
23#include <deal.II/base/point.h>
24#include <deal.II/base/table.h>
25
27#include <deal.II/grid/tria.h>
28
30
31#include <array>
32#include <limits>
33#include <map>
34
36
37#ifndef DOXYGEN
39#endif
40
64{
102 template <int dim, int spacedim>
103 void
105 const double left = 0.,
106 const double right = 1.,
107 const bool colorize = false);
108
143 template <int dim>
144 void
146 const std::vector<Point<dim>> &vertices);
147
153 template <int dim, int spacedim>
154 void
157
158
188 template <int dim, int spacedim>
189 void
191 const unsigned int repetitions,
192 const double left = 0.,
193 const double right = 1.,
194 const bool colorize = false);
195
234 template <int dim, int spacedim>
235 void
237 const Point<dim> &p1,
238 const Point<dim> &p2,
239 const bool colorize = false);
240
294 template <int dim, int spacedim>
295 void
297 const std::vector<unsigned int> &repetitions,
298 const Point<dim> &p1,
299 const Point<dim> &p2,
300 const bool colorize = false);
301
321 template <int dim>
322 void
324 const std::vector<std::vector<double>> &step_sizes,
325 const Point<dim> &p_1,
326 const Point<dim> &p_2,
327 const bool colorize = false);
328
343 template <int dim>
344 void
346 const std::vector<std::vector<double>> &spacing,
347 const Point<dim> &p,
348 const Table<dim, types::material_id> &material_id,
349 const bool colorize = false);
350
378 template <int dim, int spacedim>
379 void
381 const std::vector<unsigned int> &holes);
382
432 template <int dim>
433 void
435 const double inner_radius = 0.4,
436 const double outer_radius = 1.,
437 const double pad_bottom = 2.,
438 const double pad_top = 2.,
439 const double pad_left = 1.,
440 const double pad_right = 1.,
441 const Point<dim> &center = Point<dim>(),
442 const types::manifold_id polar_manifold_id = 0,
443 const types::manifold_id tfi_manifold_id = 1,
444 const double L = 1.,
445 const unsigned int n_slices = 2,
446 const bool colorize = false);
447
540 template <int dim>
541 void
543 const double shell_region_width = 0.03,
544 const unsigned int n_shells = 2,
545 const double skewness = 2.0,
546 const bool colorize = false);
547
571 template <int dim, int spacedim>
572 void
574 const std::vector<Point<spacedim>> &vertices,
575 const bool colorize = false);
576
596 template <int dim>
597 void
599 const Point<dim> (&corners)[dim],
600 const bool colorize = false);
601
620 template <int dim>
621 void
623 const Point<dim> (&corners)[dim],
624 const bool colorize = false);
625
641 template <int dim>
642 void
644 const unsigned int n_subdivisions,
645 const Point<dim> (&corners)[dim],
646 const bool colorize = false);
647
660 template <int dim>
661 void
663#ifndef _MSC_VER
664 const unsigned int (&n_subdivisions)[dim],
665#else
666 const unsigned int *n_subdivisions,
667#endif
668 const Point<dim> (&corners)[dim],
669 const bool colorize = false);
670
694 template <int dim, int spacedim>
695 void
697 const Point<spacedim> &origin,
698 const std::array<Tensor<1, spacedim>, dim> &edges,
699 const std::vector<unsigned int> &subdivisions = {},
700 const bool colorize = false);
701
722 template <int dim>
723 void
725 const double left = 0.,
726 const double right = 1.,
727 const double thickness = 1.,
728 const bool colorize = false);
729
809 template <int dim>
810 void
812 const Point<dim> &center = Point<dim>(),
813 const double radius = 1.,
814 const bool attach_spherical_manifold_on_boundary_cells = false);
815
851 template <int dim>
852 void
854 const Point<dim> &center = Point<dim>(),
855 const double radius = 1.);
856
875 void
877 const unsigned int n_rotate_middle_square);
878
900 void
902 const bool face_orientation,
903 const bool face_flip,
904 const bool face_rotation,
905 const bool manipulate_left_cube);
906
907
936 template <int spacedim>
937 void
940 const double radius = 1.);
941
978 template <int dim>
979 void
981 const Point<dim> &center = Point<dim>(),
982 const double radius = 1.);
983
999 template <int dim>
1000 void
1002 const Point<dim> &center = Point<dim>(),
1003 const double radius = 1.);
1004
1031 template <int dim>
1032 void
1034 const double radius = 1.,
1035 const double half_length = 1.);
1036
1037
1075 template <int dim>
1076 void
1078 const unsigned int x_subdivisions,
1079 const double radius = 1.,
1080 const double half_length = 1.);
1081
1082
1117 template <int dim>
1118 void
1120 const double radius_0 = 1.0,
1121 const double radius_1 = 0.5,
1122 const double half_length = 1.0);
1123
1247 template <int dim, int spacedim>
1248 void
1250 const std::vector<std::pair<Point<spacedim>, double>> &openings,
1251 const std::pair<Point<spacedim>, double> &bifurcation,
1252 const double aspect_ratio = 0.5);
1253
1280 template <int dim, int spacedim>
1281 void
1283 const std::vector<unsigned int> &sizes,
1284 const bool colorize_cells = false);
1285
1321 template <int dim>
1322 void
1324 const double left = -1.,
1325 const double right = 1.,
1326 const bool colorize = false);
1327
1359 template <int dim, int spacedim>
1360 void
1362 const std::vector<unsigned int> &repetitions,
1363 const Point<dim> &bottom_left,
1364 const Point<dim> &top_right,
1365 const std::vector<int> &n_cells_to_remove);
1366
1386 template <int dim>
1387 void
1389 const double left = 0.,
1390 const double right = 1.,
1391 const bool colorize = false);
1392
1454 template <int dim, int spacedim>
1455 void
1457 const Point<spacedim> &center,
1458 const double inner_radius,
1459 const double outer_radius,
1460 const unsigned int n_cells = 0,
1461 bool colorize = false);
1462
1496 template <int dim>
1497 void
1499 const Point<dim> &inner_center,
1500 const Point<dim> &outer_center,
1501 const double inner_radius,
1502 const double outer_radius,
1503 const unsigned int n_cells);
1504
1534 template <int dim>
1535 void
1537 const Point<dim> &center,
1538 const double inner_radius,
1539 const double outer_radius,
1540 const unsigned int n_cells = 0,
1541 const bool colorize = false);
1542
1543
1569 template <int dim>
1570 void
1572 const Point<dim> &center,
1573 const double inner_radius,
1574 const double outer_radius,
1575 const unsigned int n_cells = 0,
1576 const bool colorize = false);
1577
1609 template <int dim>
1610 void
1612 const double length,
1613 const double inner_radius,
1614 const double outer_radius,
1615 const unsigned int n_radial_cells = 0,
1616 const unsigned int n_axial_cells = 0,
1617 const bool colorize = false);
1618
1670 template <int dim, int spacedim>
1671 void
1673 const double centerline_radius,
1674 const double inner_radius,
1675 const unsigned int n_cells_toroidal = 6,
1676 const double phi = 2.0 * numbers::PI);
1677
1711 template <int dim, int spacedim>
1712 void
1714 const double inner_radius = .25,
1715 const double outer_radius = .5,
1716 const double L = .5,
1717 const unsigned int repetitions = 1,
1718 const bool colorize = false);
1719
1790 template <int dim>
1791 void
1793 const Point<dim> &center,
1794 const double inner_radius = 0.125,
1795 const double outer_radius = 0.25,
1796 const unsigned int n_shells = 1,
1797 const double skewness = 0.1,
1798 const unsigned int n_cells_per_shell = 0,
1799 const bool colorize = false);
1800
1814 void
1816 const unsigned int n_cells,
1817 const unsigned int n_rotations,
1818 const double R,
1819 const double r);
1820
1856 template <int dim, int spacedim>
1857 void
1860 const std::string &grid_generator_function_name,
1861 const std::string &grid_generator_function_arguments);
1862
1913 template <int dim>
1914 void
1919 const Point<3> &interior_point = Point<3>(),
1920 const double &outer_ball_radius = 1.0);
1921
1936 void
1938 Triangulation<3> &vol_tria,
2009 template <int dim, int spacedim>
2010 void
2012 const Triangulation<dim, spacedim> &triangulation_2,
2014 const double duplicated_vertex_tolerance = 1.0e-12,
2015 const bool copy_manifold_ids = false,
2016 const bool copy_boundary_ids = false);
2017
2037 template <int dim, int spacedim>
2038 void
2040 const std::vector<const Triangulation<dim, spacedim> *> &triangulations,
2042 const double duplicated_vertex_tolerance = 1.0e-12,
2043 const bool copy_manifold_ids = false,
2044 const bool copy_boundary_ids = false);
2045
2096 template <int dim, int spacedim = dim>
2097 void
2099 const std::vector<unsigned int> &extents,
2101
2139 template <int dim, int spacedim>
2140 void
2142 const Triangulation<dim, spacedim> &triangulation_1,
2143 const Triangulation<dim, spacedim> &triangulation_2,
2145
2182 template <int dim, int spacedim>
2183 void
2185 const Triangulation<dim, spacedim> &input_triangulation,
2187 &cells_to_remove,
2189
2260 void
2262 const Triangulation<2, 2> &input,
2263 const unsigned int n_slices,
2264 const double height,
2265 Triangulation<3, 3> &result,
2266 const bool copy_manifold_ids = false,
2267 const std::vector<types::manifold_id> &manifold_priorities = {});
2268
2275 void
2277 const Triangulation<2, 2> &input,
2278 const unsigned int n_slices,
2279 const double height,
2280 Triangulation<2, 2> &result,
2281 const bool copy_manifold_ids = false,
2282 const std::vector<types::manifold_id> &manifold_priorities = {});
2283
2284
2303 void
2305 const Triangulation<2, 2> &input,
2306 const std::vector<double> &slice_coordinates,
2307 Triangulation<3, 3> &result,
2308 const bool copy_manifold_ids = false,
2309 const std::vector<types::manifold_id> &manifold_priorities = {});
2310
2317 void
2319 const Triangulation<2, 2> &input,
2320 const std::vector<double> &slice_coordinates,
2321 Triangulation<2, 2> &result,
2322 const bool copy_manifold_ids = false,
2323 const std::vector<types::manifold_id> &manifold_priorities = {});
2324
2325
2326
2366 template <int dim, int spacedim1, int spacedim2>
2367 void
2370
2410 template <int dim, int spacedim>
2411 void
2414
2415 // Doxygen will not show the function above if we include the
2416 // specialization.
2417#ifndef DOXYGEN
2421 template <int spacedim>
2422 void
2424 Triangulation<1, spacedim> &out_tria);
2425#endif
2426
2431 namespace Airfoil
2432 {
2437 struct AdditionalData
2438 {
2443 std::string airfoil_type;
2444
2452 std::string naca_id;
2453
2460
2465 double airfoil_length;
2466
2471 double height;
2472
2476 double length_b2;
2477
2493 double incline_factor;
2494
2501 double bias_factor;
2502
2506 unsigned int refinements;
2507
2511 unsigned int n_subdivision_x_0;
2512
2516 unsigned int n_subdivision_x_1;
2517
2521 unsigned int n_subdivision_x_2;
2522
2526 unsigned int n_subdivision_y;
2527
2539 unsigned int airfoil_sampling_factor;
2540
2545
2551 void
2553 };
2554
2584 template <int dim>
2585 void
2588 const AdditionalData &additional_data = AdditionalData());
2589
2590
2591
2604 template <int dim>
2605 void
2608 std::vector<GridTools::PeriodicFacePair<
2609 typename Triangulation<dim, dim>::cell_iterator>> &periodic_faces,
2610 const AdditionalData &additional_data = AdditionalData());
2611
2612 } // namespace Airfoil
2613
2628 template <int dim, int spacedim>
2629 void
2632 const std::vector<unsigned int> &repetitions,
2633 const Point<dim> &p1,
2634 const Point<dim> &p2,
2635 const bool colorize = false);
2636
2652 template <int dim, int spacedim>
2653 void
2655 const unsigned int repetitions,
2656 const double p1 = 0.0,
2657 const double p2 = 1.0,
2658 const bool colorize = false);
2659
2669#ifdef _MSC_VER
2670 // Microsoft's VC++ has a bug where it doesn't want to recognize that
2671 // an implementation (definition) of the extract_boundary_mesh function
2672 // matches a declaration. This can apparently only be avoided by
2673 // doing some contortion with the return type using the following
2674 // intermediate type. This is only used when using MS VC++ and uses
2675 // the direct way of doing it otherwise
2676 template <template <int, int> class MeshType, int dim, int spacedim>
2678 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2679 struct ExtractBoundaryMesh
2680 {
2681 using return_type =
2682 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2683 typename MeshType<dim, spacedim>::face_iterator>;
2684 };
2685#endif
2686
2784 template <template <int, int> class MeshType, int dim, int spacedim>
2786 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2787#ifdef DOXYGEN
2788 return_type
2789#else
2790# ifndef _MSC_VER
2791 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2792 typename MeshType<dim, spacedim>::face_iterator>
2793# else
2794 typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
2795# endif
2796#endif
2797 extract_boundary_mesh(const MeshType<dim, spacedim> &volume_mesh,
2798 MeshType<dim - 1, spacedim> &surface_mesh,
2799 const std::set<types::boundary_id> &boundary_ids =
2800 std::set<types::boundary_id>());
2801
2819 int,
2820 << "The number of repetitions " << arg1 << " must be >=1.");
2825 int,
2826 << "The vector of repetitions must have " << arg1
2827 << " elements.");
2828
2833 "The input to this function is oriented in a way that will"
2834 " cause all cells to have negative measure.");
2837#ifndef DOXYGEN
2838 // These functions are only implemented with specializations; declare them
2839 // here
2840 template <>
2841 void
2843 const double,
2844 const double,
2845 const double,
2846 const unsigned int,
2847 const bool);
2848
2849 template <>
2850 void
2852 const double,
2853 const double,
2854 const double,
2855 const unsigned int,
2856 const bool);
2857
2858 template <>
2859 void
2861 const double,
2862 const double,
2863 const double,
2864 const unsigned int,
2865 const bool);
2866
2867 template <>
2868 void
2870 const double,
2871 const unsigned int,
2872 const double,
2873 const bool);
2874
2875 template <>
2876 void
2878 const double,
2879 const unsigned int,
2880 const double,
2881 const bool);
2882
2883 template <>
2884 void
2886 const double,
2887 const unsigned int,
2888 const double,
2889 const bool);
2890
2891
2892
2893#endif
2894} // namespace GridGenerator
2895
2896
2897
2899
2900#endif
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
Point< 3 > center
Point< 3 > vertices[4]
bool colorize
Definition grid_out.cc:4625
#define DeclException0(Exception0)
Definition exceptions.h:471
static ::ExceptionBase & ExcInvalidInputOrientation()
static ::ExceptionBase & ExcInvalidRepetitions(int arg1)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
static ::ExceptionBase & ExcInvalidRadii()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516
static ::ExceptionBase & ExcInvalidRepetitionsDimension(int arg1)
void create_triangulation(Triangulation< dim, dim > &tria, const AdditionalData &additional_data=AdditionalData())
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 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.
void hyper_cube_with_cylindrical_hole(Triangulation< dim, spacedim > &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 hyper_ball_balanced(Triangulation< dim > &tria, const Point< dim > &center=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 > &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 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_shell(Triangulation< dim, spacedim > &tria, const Point< spacedim > &center, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, 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 > &center=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 > &center, 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 > &center=Point< dim >(), const double radius=1.)
void cheese(Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &holes)
Rectangular domain with rectangular pattern of holes.
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_union_triangulation(const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result)
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 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 > &center=Point< spacedim >(), const double radius=1.)
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 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 > &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 torus(Triangulation< dim, spacedim > &tria, const double centerline_radius, const double inner_radius, const unsigned int n_cells_toroidal=6, const double phi=2.0 *numbers::PI)
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 truncated_cone(Triangulation< dim > &tria, const double radius_0=1.0, const double radius_1=0.5, const double half_length=1.0)
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, const bool copy_boundary_ids=false)
void reference_cell(Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell)
void half_hyper_ball(Triangulation< dim > &tria, const Point< dim > &center=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, const bool colorize=false)
void flatten_triangulation(const Triangulation< dim, spacedim1 > &in_tria, Triangulation< dim, spacedim2 > &out_tria)
static constexpr double PI
Definition numbers.h:259
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation
void add_parameters(ParameterHandler &prm)