Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\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{
98 template <int dim, int spacedim>
99 void
101 const double left = 0.,
102 const double right = 1.,
103 const bool colorize = false);
104
132 template <int dim>
133 void
135 const std::vector<Point<dim>> &vertices);
136
142 template <int dim, int spacedim>
143 void
146
147
174 template <int dim, int spacedim>
175 void
177 const unsigned int repetitions,
178 const double left = 0.,
179 const double right = 1.,
180 const bool colorize = false);
181
217 template <int dim, int spacedim>
218 void
220 const Point<dim> &p1,
221 const Point<dim> &p2,
222 const bool colorize = false);
223
273 template <int dim, int spacedim>
274 void
276 const std::vector<unsigned int> &repetitions,
277 const Point<dim> &p1,
278 const Point<dim> &p2,
279 const bool colorize = false);
280
296 template <int dim>
297 void
299 const std::vector<std::vector<double>> &step_sizes,
300 const Point<dim> &p_1,
301 const Point<dim> &p_2,
302 const bool colorize = false);
303
314 template <int dim>
315 void
317 const std::vector<std::vector<double>> &spacing,
318 const Point<dim> &p,
319 const Table<dim, types::material_id> &material_id,
320 const bool colorize = false);
321
346 template <int dim, int spacedim>
347 void
349 const std::vector<unsigned int> &holes);
350
394 template <int dim>
395 void
397 const double inner_radius = 0.4,
398 const double outer_radius = 1.,
399 const double pad_bottom = 2.,
400 const double pad_top = 2.,
401 const double pad_left = 1.,
402 const double pad_right = 1.,
403 const Point<dim> &center = Point<dim>(),
404 const types::manifold_id polar_manifold_id = 0,
405 const types::manifold_id tfi_manifold_id = 1,
406 const double L = 1.,
407 const unsigned int n_slices = 2,
408 const bool colorize = false);
409
498 template <int dim>
499 void
501 const double shell_region_width = 0.03,
502 const unsigned int n_shells = 2,
503 const double skewness = 2.0,
504 const bool colorize = false);
505
526 template <int dim, int spacedim>
527 void
529 const std::vector<Point<spacedim>> &vertices,
530 const bool colorize = false);
531
548 template <int dim>
549 void
551 const Point<dim> (&corners)[dim],
552 const bool colorize = false);
553
569 template <int dim>
570 void
572 const Point<dim> (&corners)[dim],
573 const bool colorize = false);
574
586 template <int dim>
587 void
589 const unsigned int n_subdivisions,
590 const Point<dim> (&corners)[dim],
591 const bool colorize = false);
592
601 template <int dim>
602 void
604#ifndef _MSC_VER
605 const unsigned int (&n_subdivisions)[dim],
606#else
607 const unsigned int *n_subdivisions,
608#endif
609 const Point<dim> (&corners)[dim],
610 const bool colorize = false);
611
635 template <int dim, int spacedim>
636 void
638 const Point<spacedim> &origin,
639 const std::array<Tensor<1, spacedim>, dim> &edges,
640 const std::vector<unsigned int> &subdivisions = {},
641 const bool colorize = false);
642
659 template <int dim>
660 void
662 const double left = 0.,
663 const double right = 1.,
664 const double thickness = 1.,
665 const bool colorize = false);
666
743 template <int dim>
744 void
746 const Point<dim> &center = Point<dim>(),
747 const double radius = 1.,
748 const bool attach_spherical_manifold_on_boundary_cells = false);
749
782 template <int dim>
783 void
785 const Point<dim> &center = Point<dim>(),
786 const double radius = 1.);
787
806 void
808 const unsigned int n_rotate_middle_square);
809
831 void
833 const bool face_orientation,
834 const bool face_flip,
835 const bool face_rotation,
836 const bool manipulate_left_cube);
837
838
864 template <int spacedim>
865 void
868 const double radius = 1.);
869
903 template <int dim>
904 void
906 const Point<dim> &center = Point<dim>(),
907 const double radius = 1.);
908
921 template <int dim>
922 void
924 const Point<dim> &center = Point<dim>(),
925 const double radius = 1.);
926
950 template <int dim>
951 void
953 const double radius = 1.,
954 const double half_length = 1.);
955
956
990 template <int dim>
991 void
993 const unsigned int x_subdivisions,
994 const double radius = 1.,
995 const double half_length = 1.);
996
997
1028 template <int dim>
1029 void
1031 const double radius_0 = 1.0,
1032 const double radius_1 = 0.5,
1033 const double half_length = 1.0);
1034
1158 template <int dim, int spacedim>
1159 void
1161 const std::vector<std::pair<Point<spacedim>, double>> &openings,
1162 const std::pair<Point<spacedim>, double> &bifurcation,
1163 const double aspect_ratio = 0.5);
1164
1191 template <int dim, int spacedim>
1192 void
1194 const std::vector<unsigned int> &sizes,
1195 const bool colorize_cells = false);
1196
1232 template <int dim>
1233 void
1235 const double left = -1.,
1236 const double right = 1.,
1237 const bool colorize = false);
1238
1270 template <int dim, int spacedim>
1271 void
1273 const std::vector<unsigned int> &repetitions,
1274 const Point<dim> &bottom_left,
1275 const Point<dim> &top_right,
1276 const std::vector<int> &n_cells_to_remove);
1277
1297 template <int dim>
1298 void
1300 const double left = 0.,
1301 const double right = 1.,
1302 const bool colorize = false);
1303
1361 template <int dim>
1362 void
1364 const Point<dim> &center,
1365 const double inner_radius,
1366 const double outer_radius,
1367 const unsigned int n_cells = 0,
1368 bool colorize = false);
1369
1403 template <int dim>
1404 void
1406 const Point<dim> &inner_center,
1407 const Point<dim> &outer_center,
1408 const double inner_radius,
1409 const double outer_radius,
1410 const unsigned int n_cells);
1411
1441 template <int dim>
1442 void
1444 const Point<dim> &center,
1445 const double inner_radius,
1446 const double outer_radius,
1447 const unsigned int n_cells = 0,
1448 const bool colorize = false);
1449
1450
1476 template <int dim>
1477 void
1479 const Point<dim> &center,
1480 const double inner_radius,
1481 const double outer_radius,
1482 const unsigned int n_cells = 0,
1483 const bool colorize = false);
1484
1516 template <int dim>
1517 void
1519 const double length,
1520 const double inner_radius,
1521 const double outer_radius,
1522 const unsigned int n_radial_cells = 0,
1523 const unsigned int n_axial_cells = 0,
1524 const bool colorize = false);
1525
1575 template <int dim, int spacedim>
1576 void
1578 const double R,
1579 const double r,
1580 const unsigned int n_cells_toroidal = 6,
1581 const double phi = 2.0 * numbers::PI);
1582
1612 template <int dim>
1613 void
1615 const double inner_radius = .25,
1616 const double outer_radius = .5,
1617 const double L = .5,
1618 const unsigned int repetitions = 1,
1619 const bool colorize = false);
1620
1691 template <int dim>
1692 void
1694 const Point<dim> &center,
1695 const double inner_radius = 0.125,
1696 const double outer_radius = 0.25,
1697 const unsigned int n_shells = 1,
1698 const double skewness = 0.1,
1699 const unsigned int n_cells_per_shell = 0,
1700 const bool colorize = false);
1701
1715 void
1717 const unsigned int n_cells,
1718 const unsigned int n_rotations,
1719 const double R,
1720 const double r);
1721
1757 template <int dim, int spacedim>
1758 void
1761 const std::string &grid_generator_function_name,
1762 const std::string &grid_generator_function_arguments);
1763
1814 template <int dim>
1815 void
1820 const Point<3> &interior_point = Point<3>(),
1821 const double &outer_ball_radius = 1.0);
1822
1837 void
1839 Triangulation<3> &vol_tria,
1910 template <int dim, int spacedim>
1911 void
1913 const Triangulation<dim, spacedim> &triangulation_2,
1915 const double duplicated_vertex_tolerance = 1.0e-12,
1916 const bool copy_manifold_ids = false,
1917 const bool copy_boundary_ids = false);
1918
1934 template <int dim, int spacedim>
1935 void
1937 const std::vector<const Triangulation<dim, spacedim> *> &triangulations,
1939 const double duplicated_vertex_tolerance = 1.0e-12,
1940 const bool copy_manifold_ids = false,
1941 const bool copy_boundary_ids = false);
1942
1990 template <int dim, int spacedim = dim>
1991 void
1993 const std::vector<unsigned int> &extents,
1995
2033 template <int dim, int spacedim>
2034 void
2036 const Triangulation<dim, spacedim> &triangulation_1,
2037 const Triangulation<dim, spacedim> &triangulation_2,
2039
2076 template <int dim, int spacedim>
2077 void
2079 const Triangulation<dim, spacedim> &input_triangulation,
2081 &cells_to_remove,
2083
2151 void
2153 const Triangulation<2, 2> &input,
2154 const unsigned int n_slices,
2155 const double height,
2156 Triangulation<3, 3> &result,
2157 const bool copy_manifold_ids = false,
2158 const std::vector<types::manifold_id> &manifold_priorities = {});
2159
2166 void
2168 const Triangulation<2, 2> &input,
2169 const unsigned int n_slices,
2170 const double height,
2171 Triangulation<2, 2> &result,
2172 const bool copy_manifold_ids = false,
2173 const std::vector<types::manifold_id> &manifold_priorities = {});
2174
2175
2194 void
2196 const Triangulation<2, 2> &input,
2197 const std::vector<double> &slice_coordinates,
2198 Triangulation<3, 3> &result,
2199 const bool copy_manifold_ids = false,
2200 const std::vector<types::manifold_id> &manifold_priorities = {});
2201
2208 void
2210 const Triangulation<2, 2> &input,
2211 const std::vector<double> &slice_coordinates,
2212 Triangulation<2, 2> &result,
2213 const bool copy_manifold_ids = false,
2214 const std::vector<types::manifold_id> &manifold_priorities = {});
2215
2216
2217
2254 template <int dim, int spacedim1, int spacedim2>
2255 void
2258
2293 template <int dim, int spacedim>
2294 void
2297
2298 // Doxygen will not show the function above if we include the
2299 // specialization.
2300#ifndef DOXYGEN
2304 template <int spacedim>
2305 void
2307 Triangulation<1, spacedim> &out_tria);
2308#endif
2309
2316 namespace Airfoil
2317 {
2323 {
2328 std::string airfoil_type;
2329
2337 std::string naca_id;
2338
2345
2351
2356 double height;
2357
2362
2379
2387
2391 unsigned int refinements;
2392
2396 unsigned int n_subdivision_x_0;
2397
2401 unsigned int n_subdivision_x_1;
2402
2406 unsigned int n_subdivision_x_2;
2407
2411 unsigned int n_subdivision_y;
2412
2425
2430
2436 void
2438 };
2439
2469 template <int dim>
2470 void
2473 const AdditionalData &additional_data = AdditionalData());
2474
2475
2476
2489 template <int dim>
2490 void
2493 std::vector<GridTools::PeriodicFacePair<
2494 typename Triangulation<dim, dim>::cell_iterator>> &periodic_faces,
2495 const AdditionalData &additional_data = AdditionalData());
2496
2497 } // namespace Airfoil
2498
2513 template <int dim, int spacedim>
2514 void
2517 const std::vector<unsigned int> &repetitions,
2518 const Point<dim> &p1,
2519 const Point<dim> &p2,
2520 const bool colorize = false);
2521
2537 template <int dim, int spacedim>
2538 void
2540 const unsigned int repetitions,
2541 const double p1 = 0.0,
2542 const double p2 = 1.0,
2543 const bool colorize = false);
2544
2554#ifdef _MSC_VER
2555 // Microsoft's VC++ has a bug where it doesn't want to recognize that
2556 // an implementation (definition) of the extract_boundary_mesh function
2557 // matches a declaration. This can apparently only be avoided by
2558 // doing some contortion with the return type using the following
2559 // intermediate type. This is only used when using MS VC++ and uses
2560 // the direct way of doing it otherwise
2561 template <template <int, int> class MeshType, int dim, int spacedim>
2563 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2564 struct ExtractBoundaryMesh
2565 {
2566 using return_type =
2567 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2568 typename MeshType<dim, spacedim>::face_iterator>;
2569 };
2570#endif
2571
2669 template <template <int, int> class MeshType, int dim, int spacedim>
2671 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2672#ifdef DOXYGEN
2673 return_type
2674#else
2675# ifndef _MSC_VER
2676 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2677 typename MeshType<dim, spacedim>::face_iterator>
2678# else
2679 typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
2680# endif
2681#endif
2682 extract_boundary_mesh(const MeshType<dim, spacedim> &volume_mesh,
2683 MeshType<dim - 1, spacedim> &surface_mesh,
2684 const std::set<types::boundary_id> &boundary_ids =
2685 std::set<types::boundary_id>());
2686
2704 int,
2705 << "The number of repetitions " << arg1 << " must be >=1.");
2710 int,
2711 << "The vector of repetitions must have " << arg1
2712 << " elements.");
2713
2718 "The input to this function is oriented in a way that will"
2719 " cause all cells to have negative measure.");
2722#ifndef DOXYGEN
2723 // These functions are only implemented with specializations; declare them
2724 // here
2725 template <>
2726 void
2728 const double,
2729 const double,
2730 const double,
2731 const unsigned int,
2732 const bool);
2733
2734 template <>
2735 void
2737 const double,
2738 const double,
2739 const double,
2740 const unsigned int,
2741 const bool);
2742
2743 template <>
2744 void
2746 const double,
2747 const double,
2748 const double,
2749 const unsigned int,
2750 const bool);
2751
2752 template <>
2753 void
2755 const double,
2756 const unsigned int,
2757 const double,
2758 const bool);
2759
2760 template <>
2761 void
2763 const double,
2764 const unsigned int,
2765 const double,
2766 const bool);
2767
2768 template <>
2769 void
2771 const double,
2772 const unsigned int,
2773 const double,
2774 const bool);
2775
2776
2777
2778#endif
2779} // namespace GridGenerator
2780
2781
2782
2784
2785#endif
Definition point.h:111
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
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_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_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 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 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)
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.)
spacedim MeshType< dim - 1, spacedim > const std::set< types::boundary_id > & boundary_ids
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)
spacedim & volume_mesh
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)
spacedim MeshType< dim - 1, spacedim > & surface_mesh
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 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 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)