Reference documentation for deal.II version 9.5.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// Copyright (C) 1999 - 2023 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
24#include <deal.II/base/point.h>
25#include <deal.II/base/table.h>
26
28#include <deal.II/grid/tria.h>
29
31
32#include <array>
33#include <limits>
34#include <map>
35
37
38#ifndef DOXYGEN
40#endif
41
65{
69 //** @{ */
70
99 template <int dim, int spacedim>
100 void
102 const double left = 0.,
103 const double right = 1.,
104 const bool colorize = false);
105
133 template <int dim>
134 void
136 const std::vector<Point<dim>> &vertices);
137
143 template <int dim, int spacedim>
144 void
147
148
175 template <int dim, int spacedim>
176 void
178 const unsigned int repetitions,
179 const double left = 0.,
180 const double right = 1.,
181 const bool colorize = false);
182
218 template <int dim, int spacedim>
219 void
221 const Point<dim> & p1,
222 const Point<dim> & p2,
223 const bool colorize = false);
224
274 template <int dim, int spacedim>
275 void
277 const std::vector<unsigned int> &repetitions,
278 const Point<dim> & p1,
279 const Point<dim> & p2,
280 const bool colorize = false);
281
297 template <int dim>
298 void
300 const std::vector<std::vector<double>> &step_sizes,
301 const Point<dim> & p_1,
302 const Point<dim> & p_2,
303 const bool colorize = false);
304
315 template <int dim>
316 void
318 const std::vector<std::vector<double>> &spacing,
319 const Point<dim> & p,
320 const Table<dim, types::material_id> &material_id,
321 const bool colorize = false);
322
347 template <int dim, int spacedim>
348 void
350 const std::vector<unsigned int> &holes);
351
395 template <int dim>
396 void
398 const double inner_radius = 0.4,
399 const double outer_radius = 1.,
400 const double pad_bottom = 2.,
401 const double pad_top = 2.,
402 const double pad_left = 1.,
403 const double pad_right = 1.,
404 const Point<dim> & center = Point<dim>(),
405 const types::manifold_id polar_manifold_id = 0,
406 const types::manifold_id tfi_manifold_id = 1,
407 const double L = 1.,
408 const unsigned int n_slices = 2,
409 const bool colorize = false);
410
499 template <int dim>
500 void
502 const double shell_region_width = 0.03,
503 const unsigned int n_shells = 2,
504 const double skewness = 2.0,
505 const bool colorize = false);
506
527 template <int dim, int spacedim>
528 void
530 const std::vector<Point<spacedim>> &vertices,
531 const bool colorize = false);
532
549 template <int dim>
550 void
552 const Point<dim> (&corners)[dim],
553 const bool colorize = false);
554
570 template <int dim>
571 void
573 const Point<dim> (&corners)[dim],
574 const bool colorize = false);
575
587 template <int dim>
588 void
590 const unsigned int n_subdivisions,
591 const Point<dim> (&corners)[dim],
592 const bool colorize = false);
593
602 template <int dim>
603 void
605#ifndef _MSC_VER
606 const unsigned int (&n_subdivisions)[dim],
607#else
608 const unsigned int *n_subdivisions,
609#endif
610 const Point<dim> (&corners)[dim],
611 const bool colorize = false);
612
636 template <int dim, int spacedim>
637 void
639 const Point<spacedim> & origin,
640 const std::array<Tensor<1, spacedim>, dim> &edges,
641 const std::vector<unsigned int> &subdivisions = {},
642 const bool colorize = false);
643
660 template <int dim>
661 void
663 const double left = 0.,
664 const double right = 1.,
665 const double thickness = 1.,
666 const bool colorize = false);
667
744 template <int dim>
745 void
747 const Point<dim> & center = Point<dim>(),
748 const double radius = 1.,
749 const bool attach_spherical_manifold_on_boundary_cells = false);
750
783 template <int dim>
784 void
786 const Point<dim> & center = Point<dim>(),
787 const double radius = 1.);
788
807 void
809 const unsigned int n_rotate_middle_square);
810
832 void
834 const bool face_orientation,
835 const bool face_flip,
836 const bool face_rotation,
837 const bool manipulate_left_cube);
838
839
865 template <int spacedim>
866 void
869 const double radius = 1.);
870
904 template <int dim>
905 void
907 const Point<dim> & center = Point<dim>(),
908 const double radius = 1.);
909
922 template <int dim>
923 void
925 const Point<dim> & center = Point<dim>(),
926 const double radius = 1.);
927
951 template <int dim>
952 void
954 const double radius = 1.,
955 const double half_length = 1.);
956
957
991 template <int dim>
992 void
994 const unsigned int x_subdivisions,
995 const double radius = 1.,
996 const double half_length = 1.);
997
998
1029 template <int dim>
1030 void
1032 const double radius_0 = 1.0,
1033 const double radius_1 = 0.5,
1034 const double half_length = 1.0);
1035
1159 template <int dim, int spacedim>
1160 void
1162 const std::vector<std::pair<Point<spacedim>, double>> &openings,
1163 const std::pair<Point<spacedim>, double> &bifurcation,
1164 const double aspect_ratio = 0.5);
1165
1192 template <int dim, int spacedim>
1193 void
1195 const std::vector<unsigned int> &sizes,
1196 const bool colorize_cells = false);
1197
1233 template <int dim>
1234 void
1236 const double left = -1.,
1237 const double right = 1.,
1238 const bool colorize = false);
1239
1271 template <int dim, int spacedim>
1272 void
1274 const std::vector<unsigned int> &repetitions,
1275 const Point<dim> & bottom_left,
1276 const Point<dim> & top_right,
1277 const std::vector<int> & n_cells_to_remove);
1278
1298 template <int dim>
1299 void
1301 const double left = 0.,
1302 const double right = 1.,
1303 const bool colorize = false);
1304
1362 template <int dim>
1363 void
1365 const Point<dim> & center,
1366 const double inner_radius,
1367 const double outer_radius,
1368 const unsigned int n_cells = 0,
1369 bool colorize = false);
1370
1404 template <int dim>
1405 void
1407 const Point<dim> & inner_center,
1408 const Point<dim> & outer_center,
1409 const double inner_radius,
1410 const double outer_radius,
1411 const unsigned int n_cells);
1412
1442 template <int dim>
1443 void
1445 const Point<dim> & center,
1446 const double inner_radius,
1447 const double outer_radius,
1448 const unsigned int n_cells = 0,
1449 const bool colorize = false);
1450
1451
1477 template <int dim>
1478 void
1480 const Point<dim> & center,
1481 const double inner_radius,
1482 const double outer_radius,
1483 const unsigned int n_cells = 0,
1484 const bool colorize = false);
1485
1512 template <int dim>
1513 void
1515 const double length,
1516 const double inner_radius,
1517 const double outer_radius,
1518 const unsigned int n_radial_cells = 0,
1519 const unsigned int n_axial_cells = 0);
1520
1570 template <int dim, int spacedim>
1571 void
1573 const double R,
1574 const double r,
1575 const unsigned int n_cells_toroidal = 6,
1576 const double phi = 2.0 * numbers::PI);
1577
1607 template <int dim>
1608 void
1610 const double inner_radius = .25,
1611 const double outer_radius = .5,
1612 const double L = .5,
1613 const unsigned int repetitions = 1,
1614 const bool colorize = false);
1615
1686 template <int dim>
1687 void
1689 const Point<dim> & center,
1690 const double inner_radius = 0.125,
1691 const double outer_radius = 0.25,
1692 const unsigned int n_shells = 1,
1693 const double skewness = 0.1,
1694 const unsigned int n_cells_per_shell = 0,
1695 const bool colorize = false);
1696
1710 void
1712 const unsigned int n_cells,
1713 const unsigned int n_rotations,
1714 const double R,
1715 const double r);
1716
1752 template <int dim, int spacedim>
1753 void
1756 const std::string & grid_generator_function_name,
1757 const std::string & grid_generator_function_arguments);
1758
1809 template <int dim>
1810 void
1815 const Point<3> &interior_point = Point<3>(),
1816 const double & outer_ball_radius = 1.0);
1817
1832 void
1834 Triangulation<3> & vol_tria,
1837 //** @} */
1838
1842 //** @{ */
1843
1905 template <int dim, int spacedim>
1906 void
1908 const Triangulation<dim, spacedim> &triangulation_2,
1910 const double duplicated_vertex_tolerance = 1.0e-12,
1911 const bool copy_manifold_ids = false,
1912 const bool copy_boundary_ids = false);
1913
1929 template <int dim, int spacedim>
1930 void
1932 const std::vector<const Triangulation<dim, spacedim> *> &triangulations,
1934 const double duplicated_vertex_tolerance = 1.0e-12,
1935 const bool copy_manifold_ids = false,
1936 const bool copy_boundary_ids = false);
1937
1985 template <int dim, int spacedim = dim>
1986 void
1988 const std::vector<unsigned int> & extents,
1990
2028 template <int dim, int spacedim>
2029 void
2031 const Triangulation<dim, spacedim> &triangulation_1,
2032 const Triangulation<dim, spacedim> &triangulation_2,
2034
2071 template <int dim, int spacedim>
2072 void
2074 const Triangulation<dim, spacedim> &input_triangulation,
2076 & cells_to_remove,
2078
2146 void
2148 const Triangulation<2, 2> & input,
2149 const unsigned int n_slices,
2150 const double height,
2151 Triangulation<3, 3> & result,
2152 const bool copy_manifold_ids = false,
2153 const std::vector<types::manifold_id> &manifold_priorities = {});
2154
2161 void
2163 const Triangulation<2, 2> & input,
2164 const unsigned int n_slices,
2165 const double height,
2166 Triangulation<2, 2> & result,
2167 const bool copy_manifold_ids = false,
2168 const std::vector<types::manifold_id> &manifold_priorities = {});
2169
2170
2189 void
2191 const Triangulation<2, 2> & input,
2192 const std::vector<double> & slice_coordinates,
2193 Triangulation<3, 3> & result,
2194 const bool copy_manifold_ids = false,
2195 const std::vector<types::manifold_id> &manifold_priorities = {});
2196
2203 void
2205 const Triangulation<2, 2> & input,
2206 const std::vector<double> & slice_coordinates,
2207 Triangulation<2, 2> & result,
2208 const bool copy_manifold_ids = false,
2209 const std::vector<types::manifold_id> &manifold_priorities = {});
2210
2211
2212
2249 template <int dim, int spacedim1, int spacedim2>
2250 void
2253
2288 template <int dim, int spacedim>
2289 void
2292
2293 // Doxygen will not show the function above if we include the
2294 // specialization.
2295#ifndef DOXYGEN
2299 template <int spacedim>
2300 void
2302 Triangulation<1, spacedim> & out_tria);
2303#endif
2304
2311 namespace Airfoil
2312 {
2318 {
2323 std::string airfoil_type;
2324
2332 std::string naca_id;
2333
2340
2346
2351 double height;
2352
2357
2374
2382
2386 unsigned int refinements;
2387
2391 unsigned int n_subdivision_x_0;
2392
2396 unsigned int n_subdivision_x_1;
2397
2401 unsigned int n_subdivision_x_2;
2402
2406 unsigned int n_subdivision_y;
2407
2420
2425
2431 void
2433 };
2434
2464 template <int dim>
2465 void
2468 const AdditionalData & additional_data = AdditionalData());
2469
2470
2471
2484 template <int dim>
2485 void
2488 std::vector<GridTools::PeriodicFacePair<
2489 typename Triangulation<dim, dim>::cell_iterator>> &periodic_faces,
2490 const AdditionalData &additional_data = AdditionalData());
2491
2492 } // namespace Airfoil
2493
2508 template <int dim, int spacedim>
2509 void
2512 const std::vector<unsigned int> &repetitions,
2513 const Point<dim> & p1,
2514 const Point<dim> & p2,
2515 const bool colorize = false);
2516
2532 template <int dim, int spacedim>
2533 void
2535 const unsigned int repetitions,
2536 const double p1 = 0.0,
2537 const double p2 = 1.0,
2538 const bool colorize = false);
2539
2540 //** @} */
2541
2547 //** @{ */
2548
2549#ifdef _MSC_VER
2550 // Microsoft's VC++ has a bug where it doesn't want to recognize that
2551 // an implementation (definition) of the extract_boundary_mesh function
2552 // matches a declaration. This can apparently only be avoided by
2553 // doing some contortion with the return type using the following
2554 // intermediate type. This is only used when using MS VC++ and uses
2555 // the direct way of doing it otherwise
2556 template <template <int, int> class MeshType, int dim, int spacedim>
2558 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2559 struct ExtractBoundaryMesh
2560 {
2561 using return_type =
2562 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2563 typename MeshType<dim, spacedim>::face_iterator>;
2564 };
2565#endif
2566
2664 template <template <int, int> class MeshType, int dim, int spacedim>
2666 (concepts::is_triangulation_or_dof_handler<MeshType<dim, spacedim>>))
2667#ifdef DOXYGEN
2668 return_type
2669#else
2670# ifndef _MSC_VER
2671 std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
2672 typename MeshType<dim, spacedim>::face_iterator>
2673# else
2674 typename ExtractBoundaryMesh<MeshType, dim, spacedim>::return_type
2675# endif
2676#endif
2677 extract_boundary_mesh(const MeshType<dim, spacedim> & volume_mesh,
2678 MeshType<dim - 1, spacedim> & surface_mesh,
2679 const std::set<types::boundary_id> &boundary_ids =
2680 std::set<types::boundary_id>());
2681
2682 //** @} */
2683
2684
2688 //** @{ */
2689
2690
2699 int,
2700 << "The number of repetitions " << arg1 << " must be >=1.");
2705 int,
2706 << "The vector of repetitions must have " << arg1
2707 << " elements.");
2708
2713 "The input to this function is oriented in a way that will"
2714 " cause all cells to have negative measure.");
2715 //** @} */
2716
2717#ifndef DOXYGEN
2718 // These functions are only implemented with specializations; declare them
2719 // here
2720 template <>
2721 void
2723 const double,
2724 const double,
2725 const double,
2726 const unsigned int,
2727 const bool);
2728
2729 template <>
2730 void
2732 const double,
2733 const double,
2734 const double,
2735 const unsigned int,
2736 const bool);
2737
2738 template <>
2739 void
2741 const double,
2742 const double,
2743 const double,
2744 const unsigned int,
2745 const bool);
2746
2747 template <>
2748 void
2750 const double,
2751 const unsigned int,
2752 const double,
2753 const bool);
2754
2755 template <>
2756 void
2758 const double,
2759 const unsigned int,
2760 const double,
2761 const bool);
2762
2763 template <>
2764 void
2766 const double,
2767 const unsigned int,
2768 const double,
2769 const bool);
2770
2771
2772
2773#endif
2774} // namespace GridGenerator
2775
2776
2777
2779
2780#endif
Definition point.h:112
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:472
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:160
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:473
Point< 3 > center
Point< 3 > vertices[4]
bool colorize
Definition grid_out.cc:4617
#define DeclException0(Exception0)
Definition exceptions.h:465
static ::ExceptionBase & ExcInvalidInputOrientation()
static ::ExceptionBase & ExcInvalidRepetitions(int arg1)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:488
static ::ExceptionBase & ExcInvalidRadii()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:510
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 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 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)
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 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 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)
const ::Triangulation< dim, spacedim > & tria