deal.II version 9.7.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 - 2025 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
29#include <deal.II/grid/tria.h>
30
31#include <array>
32#include <limits>
33#include <map>
34#include <set>
35
36
38
39#ifndef DOXYGEN
41#endif
42
66{
71
104 template <int dim, int spacedim>
105 void
107 const double left = 0.,
108 const double right = 1.,
109 const bool colorize = false);
110
145 template <int dim>
146 void
148 const std::vector<Point<dim>> &vertices);
149
155 template <int dim, int spacedim>
156 void
159
160
190 template <int dim, int spacedim>
191 void
193 const unsigned int repetitions,
194 const double left = 0.,
195 const double right = 1.,
196 const bool colorize = false);
197
236 template <int dim, int spacedim>
237 void
239 const Point<dim> &p1,
240 const Point<dim> &p2,
241 const bool colorize = false);
242
296 template <int dim, int spacedim>
297 void
299 const std::vector<unsigned int> &repetitions,
300 const Point<dim> &p1,
301 const Point<dim> &p2,
302 const bool colorize = false);
303
323 template <int dim>
324 void
326 const std::vector<std::vector<double>> &step_sizes,
327 const Point<dim> &p_1,
328 const Point<dim> &p_2,
329 const bool colorize = false);
330
345 template <int dim>
346 void
348 const std::vector<std::vector<double>> &spacing,
349 const Point<dim> &p,
350 const Table<dim, types::material_id> &material_id,
351 const bool colorize = false);
352
380 template <int dim, int spacedim>
381 void
383 const std::vector<unsigned int> &holes);
384
434 template <int dim>
435 void
437 const double inner_radius = 0.4,
438 const double outer_radius = 1.,
439 const double pad_bottom = 2.,
440 const double pad_top = 2.,
441 const double pad_left = 1.,
442 const double pad_right = 1.,
443 const Point<dim> &center = Point<dim>(),
444 const types::manifold_id polar_manifold_id = 0,
445 const types::manifold_id tfi_manifold_id = 1,
446 const double L = 1.,
447 const unsigned int n_slices = 2,
448 const bool colorize = false);
449
542 template <int dim>
543 void
545 const double shell_region_width = 0.03,
546 const unsigned int n_shells = 2,
547 const double skewness = 2.0,
548 const bool colorize = false);
549
634 template <int dim>
635 void
637 Triangulation<dim> &tria,
638 const std::vector<unsigned int> &lengths_and_heights,
639 const double depth = 1,
640 const unsigned int depth_division = 1,
641 const double shell_region_radius = 0.75,
642 const unsigned int n_shells = 2,
643 const double skewness = 2.0,
644 const bool use_transfinite_region = false,
645 const bool colorize = false);
646
670 template <int dim, int spacedim>
671 void
673 const std::vector<Point<spacedim>> &vertices,
674 const bool colorize = false);
675
695 template <int dim>
696 void
698 const Point<dim> (&corners)[dim],
699 const bool colorize = false);
700
719 template <int dim>
720 void
722 const Point<dim> (&corners)[dim],
723 const bool colorize = false);
724
740 template <int dim>
741 void
743 const unsigned int n_subdivisions,
744 const Point<dim> (&corners)[dim],
745 const bool colorize = false);
746
759 template <int dim>
760 void
762#ifndef _MSC_VER
763 const unsigned int (&n_subdivisions)[dim],
764#else
765 const unsigned int *n_subdivisions,
766#endif
767 const Point<dim> (&corners)[dim],
768 const bool colorize = false);
769
793 template <int dim, int spacedim>
794 void
796 const Point<spacedim> &origin,
797 const std::array<Tensor<1, spacedim>, dim> &edges,
798 const std::vector<unsigned int> &subdivisions = {},
799 const bool colorize = false);
800
821 template <int dim>
822 void
824 const double left = 0.,
825 const double right = 1.,
826 const double thickness = 1.,
827 const bool colorize = false);
828
911 template <int dim, int spacedim>
912 void
914 const Point<spacedim> &center = {},
915 const double radius = 1.,
916 const bool attach_spherical_manifold_on_boundary_cells = false);
917
953 template <int dim>
954 void
956 const Point<dim> &center = Point<dim>(),
957 const double radius = 1.);
958
977 void
979 const unsigned int n_rotate_middle_square);
980
1002 void
1004 const bool face_orientation,
1005 const bool face_flip,
1006 const bool face_rotation,
1007 const bool manipulate_left_cube);
1008
1009
1038 template <int spacedim>
1039 void
1041 const Point<spacedim> &center = Point<spacedim>(),
1042 const double radius = 1.);
1043
1080 template <int dim>
1081 void
1083 const Point<dim> &center = Point<dim>(),
1084 const double radius = 1.);
1085
1101 template <int dim>
1102 void
1104 const Point<dim> &center = Point<dim>(),
1105 const double radius = 1.);
1106
1133 template <int dim>
1134 void
1136 const double radius = 1.,
1137 const double half_length = 1.);
1138
1139
1177 template <int dim>
1178 void
1180 const unsigned int x_subdivisions,
1181 const double radius = 1.,
1182 const double half_length = 1.);
1183
1184
1219 template <int dim>
1220 void
1222 const double radius_0 = 1.0,
1223 const double radius_1 = 0.5,
1224 const double half_length = 1.0);
1225
1349 template <int dim, int spacedim>
1350 void
1352 const std::vector<std::pair<Point<spacedim>, double>> &openings,
1353 const std::pair<Point<spacedim>, double> &bifurcation,
1354 const double aspect_ratio = 0.5);
1355
1382 template <int dim, int spacedim>
1383 void
1385 const std::vector<unsigned int> &sizes,
1386 const bool colorize_cells = false);
1387
1423 template <int dim>
1424 void
1426 const double left = -1.,
1427 const double right = 1.,
1428 const bool colorize = false);
1429
1461 template <int dim, int spacedim>
1462 void
1464 const std::vector<unsigned int> &repetitions,
1465 const Point<dim> &bottom_left,
1466 const Point<dim> &top_right,
1467 const std::vector<int> &n_cells_to_remove);
1468
1488 template <int dim>
1489 void
1491 const double left = 0.,
1492 const double right = 1.,
1493 const bool colorize = false);
1494
1556 template <int dim, int spacedim>
1557 void
1559 const Point<spacedim> &center,
1560 const double inner_radius,
1561 const double outer_radius,
1562 const unsigned int n_cells = 0,
1563 bool colorize = false);
1564
1598 template <int dim>
1599 void
1601 const Point<dim> &inner_center,
1602 const Point<dim> &outer_center,
1603 const double inner_radius,
1604 const double outer_radius,
1605 const unsigned int n_cells);
1606
1636 template <int dim>
1637 void
1639 const Point<dim> &center,
1640 const double inner_radius,
1641 const double outer_radius,
1642 const unsigned int n_cells = 0,
1643 const bool colorize = false);
1644
1645
1671 template <int dim>
1672 void
1674 const Point<dim> &center,
1675 const double inner_radius,
1676 const double outer_radius,
1677 const unsigned int n_cells = 0,
1678 const bool colorize = false);
1679
1711 template <int dim>
1712 void
1714 const double length,
1715 const double inner_radius,
1716 const double outer_radius,
1717 const unsigned int n_radial_cells = 0,
1718 const unsigned int n_axial_cells = 0,
1719 const bool colorize = false);
1720
1772 template <int dim, int spacedim>
1773 void
1775 const double centerline_radius,
1776 const double inner_radius,
1777 const unsigned int n_cells_toroidal = 6,
1778 const double phi = 2.0 * numbers::PI);
1779
1813 template <int dim, int spacedim>
1814 void
1816 const double inner_radius = .25,
1817 const double outer_radius = .5,
1818 const double L = .5,
1819 const unsigned int repetitions = 1,
1820 const bool colorize = false);
1821
1892 template <int dim>
1893 void
1895 const Point<dim> &center,
1896 const double inner_radius = 0.125,
1897 const double outer_radius = 0.25,
1898 const unsigned int n_shells = 1,
1899 const double skewness = 0.1,
1900 const unsigned int n_cells_per_shell = 0,
1901 const bool colorize = false);
1902
1916 void
1918 const unsigned int n_cells,
1919 const unsigned int n_rotations,
1920 const double R,
1921 const double r);
1922
1958 template <int dim, int spacedim>
1959 void
1962 const std::string &grid_generator_function_name,
1963 const std::string &grid_generator_function_arguments);
1964
2015 template <int dim>
2016 void
2021 const Point<3> &interior_point = Point<3>(),
2022 const double &outer_ball_radius = 1.0);
2023
2038 void
2040 Triangulation<3> &vol_tria,
2044
2049
2111 template <int dim, int spacedim>
2112 void
2114 const Triangulation<dim, spacedim> &triangulation_2,
2116 const double duplicated_vertex_tolerance = 1.0e-12,
2117 const bool copy_manifold_ids = false,
2118 const bool copy_boundary_ids = false);
2119
2139 template <int dim, int spacedim>
2140 void
2142 const std::vector<const Triangulation<dim, spacedim> *> &triangulations,
2144 const double duplicated_vertex_tolerance = 1.0e-12,
2145 const bool copy_manifold_ids = false,
2146 const bool copy_boundary_ids = false);
2147
2198 template <int dim, int spacedim = dim>
2199 void
2201 const std::vector<unsigned int> &extents,
2203
2241 template <int dim, int spacedim>
2242 void
2244 const Triangulation<dim, spacedim> &triangulation_1,
2245 const Triangulation<dim, spacedim> &triangulation_2,
2247
2284 template <int dim, int spacedim>
2285 void
2287 const Triangulation<dim, spacedim> &input_triangulation,
2289 &cells_to_remove,
2291
2362 void
2364 const Triangulation<2, 2> &input,
2365 const unsigned int n_slices,
2366 const double height,
2367 Triangulation<3, 3> &result,
2368 const bool copy_manifold_ids = false,
2369 const std::vector<types::manifold_id> &manifold_priorities = {});
2370
2377 void
2379 const Triangulation<2, 2> &input,
2380 const unsigned int n_slices,
2381 const double height,
2382 Triangulation<2, 2> &result,
2383 const bool copy_manifold_ids = false,
2384 const std::vector<types::manifold_id> &manifold_priorities = {});
2385
2386
2405 void
2407 const Triangulation<2, 2> &input,
2408 const std::vector<double> &slice_coordinates,
2409 Triangulation<3, 3> &result,
2410 const bool copy_manifold_ids = false,
2411 const std::vector<types::manifold_id> &manifold_priorities = {});
2412
2419 void
2421 const Triangulation<2, 2> &input,
2422 const std::vector<double> &slice_coordinates,
2423 Triangulation<2, 2> &result,
2424 const bool copy_manifold_ids = false,
2425 const std::vector<types::manifold_id> &manifold_priorities = {});
2426
2427
2428
2468 template <int dim, int spacedim1, int spacedim2>
2469 void
2472
2526 template <int dim, int spacedim>
2527 void
2530 const unsigned int n_divisions = (dim == 2 ?
2531 8u :
2532 24u));
2533
2556 template <int dim, int spacedim>
2557 void
2560
2565 namespace Airfoil
2566 {
2571 struct AdditionalData
2572 {
2577 std::string airfoil_type;
2578
2586 std::string naca_id;
2587
2594
2599 double airfoil_length;
2600
2605 double height;
2606
2610 double length_b2;
2611
2627 double incline_factor;
2628
2635 double bias_factor;
2636
2640 unsigned int refinements;
2641
2645 unsigned int n_subdivision_x_0;
2646
2650 unsigned int n_subdivision_x_1;
2651
2655 unsigned int n_subdivision_x_2;
2656
2660 unsigned int n_subdivision_y;
2661
2673 unsigned int airfoil_sampling_factor;
2674
2679
2685 void
2687 };
2688
2718 template <int dim>
2719 void
2722 const AdditionalData &additional_data = AdditionalData());
2723
2724
2725
2738 template <int dim>
2739 void
2742 std::vector<GridTools::PeriodicFacePair<
2743 typename Triangulation<dim, dim>::cell_iterator>> &periodic_faces,
2744 const AdditionalData &additional_data = AdditionalData());
2745
2746 } // namespace Airfoil
2747
2764 template <int dim, int spacedim>
2765 void
2768 const std::vector<unsigned int> &repetitions,
2769 const Point<dim> &p1,
2770 const Point<dim> &p2,
2771 const bool colorize = false);
2772
2790 template <int dim, int spacedim>
2791 void
2793 const unsigned int repetitions,
2794 const double p1 = 0.0,
2795 const double p2 = 1.0,
2796 const bool colorize = false);
2797
2799
2806
2807#ifdef _MSC_VER
2808 // Microsoft's VC++ has a bug where it doesn't want to recognize that
2809 // an implementation (definition) of the extract_boundary_mesh function
2810 // matches a declaration. This can apparently only be avoided by
2811 // doing some contortion with the return type using the following
2812 // intermediate type. This is only used when using MS VC++ and uses
2813 // the direct way of doing it otherwise
2814 template <template <int, int> class MeshType, int dim, int spacedim>
2816 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2817 struct ExtractBoundaryMesh
2818 {
2819 using return_type =
2820 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2821 typename MeshType<dim, spacedim>::face_iterator>;
2822 };
2823#endif
2824
2922 template <template <int, int> class MeshType, int dim, int spacedim>
2924 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2925#ifdef DOXYGEN
2926 return_type
2927#else
2928# ifndef _MSC_VER
2929 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2930 typename MeshType<dim, spacedim>::face_iterator>
2931# else
2932 typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
2933# endif
2934#endif
2935 extract_boundary_mesh(const MeshType<dim, spacedim> &volume_mesh,
2936 MeshType<dim - 1, spacedim> &surface_mesh,
2937 const std::set<types::boundary_id> &boundary_ids =
2938 std::set<types::boundary_id>());
2939
2941
2942
2947
2948
2957 int,
2958 << "The number of repetitions " << arg1 << " must be >=1.");
2963 int,
2964 << "The vector of repetitions must have " << arg1
2965 << " elements.");
2966
2971 "The input to this function is oriented in a way that will"
2972 " cause all cells to have negative measure.");
2974
2975#ifndef DOXYGEN
2976 // These functions are only implemented with specializations; declare them
2977 // here
2978 template <>
2979 void
2981 const double,
2982 const double,
2983 const double,
2984 const unsigned int,
2985 const bool);
2986
2987 template <>
2988 void
2990 const double,
2991 const double,
2992 const double,
2993 const unsigned int,
2994 const bool);
2995
2996 template <>
2997 void
2999 const double,
3000 const double,
3001 const double,
3002 const unsigned int,
3003 const bool);
3004
3005 template <>
3006 void
3008 const double,
3009 const unsigned int,
3010 const double,
3011 const bool);
3012
3013 template <>
3014 void
3016 const double,
3017 const unsigned int,
3018 const double,
3019 const bool);
3020
3021 template <>
3022 void
3024 const double,
3025 const unsigned int,
3026 const double,
3027 const bool);
3028
3029
3030
3031#endif
3032} // namespace GridGenerator
3033
3034
3035
3037
3038#endif
Definition point.h:113
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:248
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#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)
TriaActiveIterator< CellAccessor< dim, spacedim > > active_cell_iterator
Definition tria.h:1581
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
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 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 uniform_channel_with_cylinder(Triangulation< dim > &tria, const std::vector< unsigned int > &lengths_and_heights, const double depth=1, const unsigned int depth_division=1, const double shell_region_radius=0.75, const unsigned int n_shells=2, const double skewness=2.0, const bool use_transfinite_region=false, 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 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 hyper_ball(Triangulation< dim, spacedim > &tria, const Point< spacedim > &center={}, const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=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 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 alfeld_split_of_simplex_mesh(const Triangulation< dim, spacedim > &in_tria, Triangulation< dim, spacedim > &out_tria)
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 convert_hypercube_to_simplex_mesh(const Triangulation< dim, spacedim > &in_tria, Triangulation< dim, spacedim > &out_tria, const unsigned int n_divisions=(dim==2 ? 8u :24u))
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 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)
constexpr double PI
Definition numbers.h:239
unsigned int manifold_id
Definition types.h:173
void add_parameters(ParameterHandler &prm)