deal.II version GIT relicensing-1991-gec5bc16adb 2024-10-14 08:00:00+00:00
|
Namespaces | |
namespace | Airfoil |
Functions | |
Creating meshes for basic geometries | |
template<int dim, int spacedim> | |
void | hyper_cube (Triangulation< dim, spacedim > &tria, const double left=0., const double right=1., const bool colorize=false) |
template<int dim> | |
void | simplex (Triangulation< dim, dim > &tria, const std::vector< Point< dim > > &vertices) |
template<int dim, int spacedim> | |
void | reference_cell (Triangulation< dim, spacedim > &tria, const ReferenceCell &reference_cell) |
template<int dim, int spacedim> | |
void | subdivided_hyper_cube (Triangulation< dim, spacedim > &tria, const unsigned int repetitions, const double left=0., const double right=1., const bool colorize=false) |
template<int dim, int spacedim> | |
void | hyper_rectangle (Triangulation< dim, spacedim > &tria, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false) |
template<int dim, int spacedim> | |
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) |
template<int dim> | |
void | subdivided_hyper_rectangle (Triangulation< dim > &tria, const std::vector< std::vector< double > > &step_sizes, const Point< dim > &p_1, const Point< dim > &p_2, const bool colorize=false) |
template<int dim> | |
void | subdivided_hyper_rectangle (Triangulation< dim > &tria, const std::vector< std::vector< double > > &spacing, const Point< dim > &p, const Table< dim, types::material_id > &material_id, const bool colorize=false) |
template<int dim, int spacedim> | |
void | cheese (Triangulation< dim, spacedim > &tria, const std::vector< unsigned int > &holes) |
Rectangular domain with rectangular pattern of holes. | |
template<int dim> | |
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 > ¢er=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. | |
template<int dim> | |
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) |
template<int dim, int spacedim> | |
void | general_cell (Triangulation< dim, spacedim > &tria, const std::vector< Point< spacedim > > &vertices, const bool colorize=false) |
template<int dim> | |
void | parallelogram (Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false) |
template<int dim> | |
void | parallelepiped (Triangulation< dim > &tria, const Point< dim >(&corners)[dim], const bool colorize=false) |
template<int dim> | |
void | subdivided_parallelepiped (Triangulation< dim > &tria, const unsigned int n_subdivisions, const Point< dim >(&corners)[dim], const bool colorize=false) |
template<int dim> | |
void | subdivided_parallelepiped (Triangulation< dim > &tria, const unsigned int(&n_subdivisions)[dim], const Point< dim >(&corners)[dim], const bool colorize=false) |
template<int dim, int spacedim> | |
void | subdivided_parallelepiped (Triangulation< dim, spacedim > &tria, const Point< spacedim > &origin, const std::array< Tensor< 1, spacedim >, dim > &edges, const std::vector< unsigned int > &subdivisions={}, const bool colorize=false) |
template<int dim> | |
void | enclosed_hyper_cube (Triangulation< dim > &tria, const double left=0., const double right=1., const double thickness=1., const bool colorize=false) |
template<int dim> | |
void | hyper_ball (Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1., const bool attach_spherical_manifold_on_boundary_cells=false) |
template<int dim> | |
void | hyper_ball_balanced (Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.) |
void | non_standard_orientation_mesh (Triangulation< 2 > &tria, const unsigned int n_rotate_middle_square) |
void | non_standard_orientation_mesh (Triangulation< 3 > &tria, const bool face_orientation, const bool face_flip, const bool face_rotation, const bool manipulate_left_cube) |
template<int spacedim> | |
void | hyper_sphere (Triangulation< spacedim - 1, spacedim > &tria, const Point< spacedim > ¢er=Point< spacedim >(), const double radius=1.) |
template<int dim> | |
void | quarter_hyper_ball (Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.) |
template<int dim> | |
void | half_hyper_ball (Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.) |
template<int dim> | |
void | cylinder (Triangulation< dim > &tria, const double radius=1., const double half_length=1.) |
template<int dim> | |
void | subdivided_cylinder (Triangulation< dim > &tria, const unsigned int x_subdivisions, const double radius=1., const double half_length=1.) |
template<int dim> | |
void | truncated_cone (Triangulation< dim > &tria, const double radius_0=1.0, const double radius_1=0.5, const double half_length=1.0) |
template<int dim, int spacedim> | |
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) |
template<int dim, int spacedim> | |
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. | |
template<int dim> | |
void | hyper_L (Triangulation< dim > &tria, const double left=-1., const double right=1., const bool colorize=false) |
template<int dim, int spacedim> | |
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) |
template<int dim> | |
void | hyper_cube_slit (Triangulation< dim > &tria, const double left=0., const double right=1., const bool colorize=false) |
template<int dim, int spacedim> | |
void | hyper_shell (Triangulation< dim, spacedim > &tria, const Point< spacedim > ¢er, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, bool colorize=false) |
template<int dim> | |
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) |
template<int dim> | |
void | half_hyper_shell (Triangulation< dim > &tria, const Point< dim > ¢er, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, const bool colorize=false) |
template<int dim> | |
void | quarter_hyper_shell (Triangulation< dim > &tria, const Point< dim > ¢er, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, const bool colorize=false) |
template<int dim> | |
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) |
template<int dim, int spacedim> | |
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) |
template<int dim, int spacedim> | |
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) |
template<int dim> | |
void | concentric_hyper_shells (Triangulation< dim > &triangulation, const Point< dim > ¢er, 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 | moebius (Triangulation< 3, 3 > &tria, const unsigned int n_cells, const unsigned int n_rotations, const double R, const double r) |
template<int dim, int spacedim> | |
void | generate_from_name_and_arguments (Triangulation< dim, spacedim > &tria, const std::string &grid_generator_function_name, const std::string &grid_generator_function_arguments) |
template<int dim> | |
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 | surface_mesh_to_volumetric_mesh (const Triangulation< 2, 3 > &surface_tria, Triangulation< 3 > &vol_tria, const CGALWrappers::AdditionalData< 3 > &data=CGALWrappers::AdditionalData< 3 >{}) |
Creating meshes from other meshes | |
template<int dim, int spacedim> | |
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) |
template<int dim, int spacedim> | |
void | merge_triangulations (const std::vector< const Triangulation< dim, spacedim > * > &triangulations, Triangulation< dim, spacedim > &result, const double duplicated_vertex_tolerance=1.0e-12, const bool copy_manifold_ids=false, const bool copy_boundary_ids=false) |
template<int dim, int spacedim = dim> | |
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. | |
template<int dim, int spacedim> | |
void | create_union_triangulation (const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result) |
template<int dim, int spacedim> | |
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 | 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 | extrude_triangulation (const Triangulation< 2, 2 > &input, const unsigned int n_slices, const double height, Triangulation< 2, 2 > &result, const bool copy_manifold_ids=false, const std::vector< types::manifold_id > &manifold_priorities={}) |
void | extrude_triangulation (const Triangulation< 2, 2 > &input, const std::vector< double > &slice_coordinates, Triangulation< 3, 3 > &result, const bool copy_manifold_ids=false, const std::vector< types::manifold_id > &manifold_priorities={}) |
void | extrude_triangulation (const Triangulation< 2, 2 > &input, const std::vector< double > &slice_coordinates, Triangulation< 2, 2 > &result, const bool copy_manifold_ids=false, const std::vector< types::manifold_id > &manifold_priorities={}) |
template<int dim, int spacedim1, int spacedim2> | |
void | flatten_triangulation (const Triangulation< dim, spacedim1 > &in_tria, Triangulation< dim, spacedim2 > &out_tria) |
template<int dim, int spacedim> | |
void | convert_hypercube_to_simplex_mesh (const Triangulation< dim, spacedim > &in_tria, Triangulation< dim, spacedim > &out_tria) |
template<int dim, int spacedim> | |
void | alfeld_split_of_simplex_mesh (const Triangulation< dim, spacedim > &in_tria, Triangulation< dim, spacedim > &out_tria) |
template<int dim, int spacedim> | |
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) |
template<int dim, int spacedim> | |
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) |
Creating lower-dimensional meshes | |
Created from parts of higher-dimensional meshes. | |
template<template< int, int > class MeshType, int dim, int spacedim> | |
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 >()) |
Exceptions | |
static ::ExceptionBase & | ExcInvalidRadii () |
static ::ExceptionBase & | ExcInvalidRepetitions (int arg1) |
static ::ExceptionBase & | ExcInvalidRepetitionsDimension (int arg1) |
static ::ExceptionBase & | ExcInvalidInputOrientation () |
This namespace provides a collection of functions for generating triangulations for some basic geometries.
Some of these functions receive a flag colorize
(see the glossary entry on colorization). If this is set, parts of the boundary receive different boundary indicators allowing them to be distinguished for the purpose of evaluating different boundary conditions.
If the domain is curved, each of the domain parts that should be refined by following an appropriate Manifold description will receive a different manifold indicator, and the correct Manifold descriptor will be attached to the Triangulation. Notice that if you later transform the triangulation, you have to make sure you attach the correct new Manifold to the triangulation.
void GridGenerator::hyper_cube | ( | Triangulation< dim, spacedim > & | tria, |
const double | left = 0. , |
||
const double | right = 1. , |
||
const bool | colorize = false |
||
) |
Initialize the given triangulation with a hypercube (line in 1d, square in 2d, etc) consisting of exactly one cell. The hypercube volume is the tensor product interval \([left,right]^{\text{dim}}\) in the present number of dimensions, where the limits are given as arguments. They default to zero and unity, then producing the unit hypercube.
If the argument colorize
is false, then all boundary indicators are set to zero (the default boundary indicator) for 2d and 3d. If it is true, the boundary is colorized as in hyper_rectangle(). In 1d the indicators are always colorized, see hyper_rectangle().
If dim
< spacedim
, this will create a dim
dimensional object in the first dim
coordinate directions embedded into the spacedim
dimensional space with the remaining entries set to zero. For example, a Triangulation<2,3>
will be a square in the xy plane with z=0.
See also subdivided_hyper_cube() for a coarse mesh consisting of several cells. See hyper_rectangle(), if different lengths in different ordinate directions are required.
triangulation.generate_hyper_cube(left = 0., right = 1., colorize = false)
. void GridGenerator::simplex | ( | Triangulation< dim, dim > & | tria, |
const std::vector< Point< dim > > & | vertices | ||
) |
Create a \(d\)-simplex (i.e., a triangle in 2d, or a tetrahedron in 3d) with \(d+1\) corners. Since this function dates back to a time when deal.II did not support triangular and tetrahedral cells (which it does now), the simplex created by this function and described by the input arguments is subdivided into quadrilaterals and hexahedra by adding edge, face, and simplex midpoints, resulting in a mesh that consists of \(d+1\) quadrilateral or hexahedral cells. (However, you can create a simplex geometry with a single simplex cell – i.e., a triangle or a tetrahedron – with both the GridGenerator::reference_cell() and directly with the Triangulation::create_triangulation() functions.)
The vertices
argument contains a vector with all d+1 vertices defining the corners of the simplex. They must be given in an order such that the vectors from the first vertex to each of the others form a right-handed system.
The meshes generated in two and three dimensions are:
tria | The triangulation to be created. It needs to be empty upon calling this function. |
vertices | The dim+1 corners of the simplex. |
Triangulation<2,2>
, Triangulation<3,3>
.triangulation.generate_simplex(vertices)
. void GridGenerator::reference_cell | ( | Triangulation< dim, spacedim > & | tria, |
const ReferenceCell & | reference_cell | ||
) |
Create a (coarse) grid with a single cell of the shape of the provided reference cell. This is a generalization of the hyper_cube() and simplex() functions above.
void GridGenerator::subdivided_hyper_cube | ( | Triangulation< dim, spacedim > & | tria, |
const unsigned int | repetitions, | ||
const double | left = 0. , |
||
const double | right = 1. , |
||
const bool | colorize = false |
||
) |
Same as hyper_cube(), but with the difference that not only one cell is created but each coordinate direction is subdivided into repetitions
cells. Thus, the number of cells filling the given volume is repetitionsdim
.
If dim
< spacedim
, this will create a dim
dimensional object in the first dim
coordinate directions embedded into the spacedim
dimensional space with the remaining entries set to zero. For example, a Triangulation<2,3>
will be a square in the xy plane with z=0.
tria | The triangulation to create. It needs to be empty upon calling this function. |
repetitions | An unsigned integer denoting the number of cells to generate in each direction. |
left | Lower bound for the interval used to create the hyper cube. |
right | Upper bound for the interval used to create the hyper cube. |
colorize | Assign different boundary ids if set to true. |
triangulation.generate_subdivided_hyper_cube(left = 0., right = 1.)
. void GridGenerator::hyper_rectangle | ( | Triangulation< dim, spacedim > & | tria, |
const Point< dim > & | p1, | ||
const Point< dim > & | p2, | ||
const bool | colorize = false |
||
) |
Create a coordinate-parallel brick from the two diagonally opposite corner points p1
and p2
.
If the colorize
flag is true
, then the boundary_ids
of the boundary faces are assigned, such that the lower one in x-direction
is 0, the upper one is 1. The indicators for the surfaces in y-direction
are 2 and 3, the ones for z
are 4 and 5. This corresponds to the numbers of faces of the unit square of cube as laid out in the documentation of the GeometryInfo class; see also the glossary entry on colorization. Importantly, however, in 3d colorization does not set boundary_ids
of edges, but only of faces, because each boundary edge is shared between two faces and it is not clear how the boundary id of an edge should be set in that case.
Additionally, if colorize
is true
, material ids are assigned to the cells according to the octant their center is in: being in the right half space for any coordinate direction xi adds 2i. For instance, a cell with center point (1,-1,1) yields a material id 5, assuming that the center of the hyper rectangle lies at the origin. No manifold id is set for the cells.
If dim
< spacedim
, this will create a dim
dimensional object in the first dim
coordinate directions embedded into the spacedim
dimensional space with the remaining entries set to zero. For example, a Triangulation<2,3>
will be a rectangle in the xy plane with z=0, defined by the two opposing corners p1
and p2
.
triangulation.generate_hyper_rectangle(p1, p2, colorize = false)
. void GridGenerator::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 |
||
) |
Create a coordinate-parallel brick from the two diagonally opposite corner points p1
and p2
. The number of cells in coordinate direction i
is given by the integer repetitions[i]
.
To get cells with an aspect ratio different from that of the domain, use different numbers of subdivisions, given by repetitions
, in different coordinate directions. The minimum number of subdivisions in each direction is 1.
If the colorize
flag is true
, then the boundary_ids
of the surfaces are assigned, such that the lower one in x-direction
is 0, the upper one is 1 (the left and the right vertical face). The indicators for the surfaces in y-direction
are 2 and 3, the ones for z
are 4 and 5. Additionally, material ids are assigned to the cells according to the octant their center is in: being in the right half plane for any coordinate direction xi adds 2i (see the glossary entry on colorization). For instance, the center point (1,-1,1) yields a material id 5 (this means that in 2d only material ids 0,1,2,3 are assigned independent from the number of repetitions).
Note that the colorize
flag is ignored in 1d and is assumed to always be true. That means the boundary indicator is 0 on the left and 1 on the right. See step-15 for details.
If dim
< spacedim
, this will create a dim
dimensional object in the first dim
coordinate directions embedded into the spacedim
dimensional space with the remaining entries set to zero. For example, a Triangulation<2,3>
will be a rectangle in the xy plane with z=0, defined by the two opposing corners p1
and p2
.
tria | The triangulation to be created. It needs to be empty upon calling this function. |
repetitions | A vector of dim positive values denoting the number of cells to generate in that direction. |
p1 | First corner point. |
p2 | Second corner opposite to p1 . |
colorize | Assign different boundary ids if set to true. The same comments apply as for the hyper_rectangle() function. |
triangulation.generate_subdivided_hyper_rectangle(repetitions, p1, p2, colorize = false)
. void GridGenerator::subdivided_hyper_rectangle | ( | Triangulation< dim > & | tria, |
const std::vector< std::vector< double > > & | step_sizes, | ||
const Point< dim > & | p_1, | ||
const Point< dim > & | p_2, | ||
const bool | colorize = false |
||
) |
Like the previous function. However, here the second argument does not denote the number of subdivisions in each coordinate direction, but a sequence of step sizes for each coordinate direction. The domain will therefore be subdivided into step_sizes[i].size()
cells in coordinate direction i
, with width step_sizes[i][j]
for the j
th cell.
This function is therefore the right one to generate graded meshes where cells are concentrated in certain areas, rather than a uniformly subdivided mesh as the previous function generates.
The step sizes have to add up to the dimensions of the hyper rectangle specified by the points p1
and p2
.
triangulation.generate_subdivided_steps_hyper_rectangle(step_sizes, p1, p2, colorize = false)
. void GridGenerator::subdivided_hyper_rectangle | ( | Triangulation< dim > & | tria, |
const std::vector< std::vector< double > > & | spacing, | ||
const Point< dim > & | p, | ||
const Table< dim, types::material_id > & | material_id, | ||
const bool | colorize = false |
||
) |
Like the previous function, but with the following twist: the material_id
argument is a dim-dimensional array that, for each cell, indicates which material_id should be set. In addition, and this is the major new functionality, if the material_id of a cell is (unsigned char)(-1)
, then that cell is deleted from the triangulation, i.e. the domain will have a void there.
triangulation.generate_subdivided_material_hyper_rectangle(spacing, p, material_id, colorize = false)
. void GridGenerator::cheese | ( | Triangulation< dim, spacedim > & | tria, |
const std::vector< unsigned int > & | holes | ||
) |
Rectangular domain with rectangular pattern of holes.
The domain itself is rectangular, very much as if it had been generated by subdivided_hyper_rectangle(). The argument holes
specifies how many square holes the domain should have in each coordinate direction. The total number of mesh cells in that direction is then twice this number plus one.
The number of holes in one direction must be at least one.
An example with two by three holes is
If dim
< spacedim
, this will create a dim
dimensional object in the first dim
coordinate directions embedded into the spacedim
dimensional space with the remaining entries set to zero.
tria | The triangulation to be created. It needs to be empty upon calling this function. |
holes | Positive number of holes in each of the dim directions. |
triangulation.generate_cheese(holes)
. void GridGenerator::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.
Generate a rectangular plate with an (offset) cylindrical hole. The geometry consists of 2 regions: The first is a square region with length outer_radius
and a hole of radius inner_radius
. Cells in this region will have TransfiniteInterpolationManifold with manifold id tfi_manifold_id
attached to them. Additionally, the boundary faces of the hole will be associated with a PolarManifold (in 2d) or CylindricalManifold (in 3d). The center of this region can be prescribed via center
, namely the axis of the hole will be located at center
. The second region describes the remainder of the bulk material. It is specified via padding parameters pad_bottom
, padding_top
, padding_left
and padding_right
. All cells in this region will have a FlatManifold attached to them. The final width of the plate will be padding_left + 2*outer_radius + padding_right
, while its length is padding_top + 2*outer_radius + padding_bottom
.
Here is the non-symmetric grid (after one global refinement, colored according to manifold id) in 2d and 3d, respectively:
In 3d, triangulation will be extruded in the z-direction by the total height of L
using n_slices
slices (minimum is 2).
If the colorize
flag is true
, the boundary_ids of the boundary faces are assigned such that the lower one in the x-direction is 0, and the upper one is 1 (see the glossary entry on colorization). The indicators for the surfaces in the y-direction are 2 and 3, and the ones for the z-direction are 5 and 6. The hole boundary has indicator 4.
tria
is the triangulation to be created. It needs to be empty upon calling this function.
triangulation.generate_plate_with_a_holes(inner_radius = 0.4, outer_radius = 1, pad_bottom = 2, pad_top = 2., pad_left = 1., pad_right, center = Point(), polar_manifold = 0, tfi_manifold = 1, L = 1, n_slices = 2, colorize = false)
. void GridGenerator::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 |
||
) |
Generate a grid consisting of a channel with a cylinder. This is a common benchmark for Navier-Stokes solvers. The geometry consists of a channel of size \([0, 2.2] \times [0, 0.41] \times [0, 0.41] \) (where the \(z\) dimension is omitted in 2d) with a cylinder, parallel to the \(z\) axis with diameter \(0.1\), centered at \((0.2, 0.2, 0)\). The channel has three distinct regions:
n_shells
is greater than zero, then there are that many shells centered around the cylinder, Since the cylinder is slightly offset from the center of the channel, this geometry results in vortex shedding at moderate Reynolds numbers. Here is the grid (without additional global refinement) in 2d:
and in 3d:
The resulting Triangulation uses three manifolds: a PolarManifold (in 2d) or CylindricalManifold (in 3d) with manifold id \(0\), a TransfiniteInterpolationManifold with manifold id \(1\), and a FlatManifold everywhere else. For more information on this topic see the glossary entry on manifold indicators. The cell faces on the cylinder and surrounding shells have manifold ids of \(0\), while the cell volumes adjacent to the shells (or, if they do not exist, the cylinder) have a manifold id of \(1\). Put another way: this grid uses TransfiniteInterpolationManifold to smoothly transition from the shells (generated with GridGenerator::concentric_hyper_shells) to the bulk region. All other cell volumes and faces have manifold id numbers::flat_manifold_id and use FlatManifold. All cells with id numbers::flat_manifold_id are rectangular prisms aligned with the coordinate axes.
The picture below shows part of the 2d grid (using all default arguments to this function) after two global refinements. The cells with manifold id \(0\) are orange (the polar manifold id), cells with manifold id \(1\) are yellow (the transfinite interpolation manifold id), and the ones with manifold id numbers::flat_manifold_id are cyan:
tria | Triangulation to be created. Must be empty upon calling this function. |
shell_region_width | Width of the layer of shells around the cylinder. This value should be between \(0\) and \(0.05\); the default value is \(0.03\). |
n_shells | Number of shells to use in the shell layer. |
skewness | Parameter controlling how close the shells are to the cylinder: see the mathematical definition given in GridGenerator::concentric_hyper_shells. |
colorize | If true , then assign different boundary ids to different parts of the boundary. For more information on boundary indicators see this glossary entry. The left boundary (at \(x = 0\)) is assigned an id of \(0\), the right boundary (at \(x = 2.2\)) is assigned an id of \(1\); the boundary of the obstacle in the middle (i.e., the circle in 2d or the cylinder walls in 3d) is assigned an id of \(2\), and the channel walls are assigned an id of \(3\). |
See the original paper for more information:
triangulation.generate_channel_with_cylinder(shell_region_width = 0.03, n_shellsa = 2, skewness = 2., colorize = false)
. void GridGenerator::general_cell | ( | Triangulation< dim, spacedim > & | tria, |
const std::vector< Point< spacedim > > & | vertices, | ||
const bool | colorize = false |
||
) |
A general dim
-dimensional cell (a segment if dim is 1, a quadrilateral if dim
is 2, or a hexahedron if dim
is 3) immersed in a spacedim
-dimensional space. It is the responsibility of the user to provide the vertices in the right order (see the documentation of the GeometryInfo class) because the vertices are stored in the same order as they are given. It is also important to make sure that the volume of the cell is positive.
If the argument colorize
is false, then all boundary indicators are set to zero for 2d and 3d. If it is true, the boundary is colorized as in hyper_rectangle() (see the glossary entry on colorization). In 1d the indicators are always colorized, see hyper_rectangle().
tria | The triangulation that will be created |
vertices | The 2^dim vertices of the cell |
colorize | If true, set different boundary ids. |
triangulation.generate_general_cell(vertices, colorize = false)
. void GridGenerator::parallelogram | ( | Triangulation< dim > & | tria, |
const Point< dim >(&) | corners[dim], | ||
const bool | colorize = false |
||
) |
A parallelogram. The first corner point is the origin. The next dim
vertices are the ones given in the second argument and the last vertex will be the sum of the two vectors connecting the origin to those points. Colorizing is done in the same way as in hyper_rectangle().
tria | The triangulation to be created. It needs to be empty upon calling this function. |
corners | Second and third vertices of the parallelogram. |
colorize | Assign different boundary ids if true. (see the glossary entry on colorization). |
triangulation.generate_parallelogram(vertices, colorize = false)
. void GridGenerator::parallelepiped | ( | Triangulation< dim > & | tria, |
const Point< dim >(&) | corners[dim], | ||
const bool | colorize = false |
||
) |
A parallelepiped. The first corner point is the origin. The dim
adjacent points are vectors describing the edges of the parallelepiped with respect to the origin. Additional points are sums of these dim vectors. Colorizing is done according to hyper_rectangle().
corners
will no longer refer to the same triangulation.triangulation.generate_parallelepiped(vertices, colorize = false)
. void GridGenerator::subdivided_parallelepiped | ( | Triangulation< dim > & | tria, |
const unsigned int | n_subdivisions, | ||
const Point< dim >(&) | corners[dim], | ||
const bool | colorize = false |
||
) |
A subdivided parallelepiped. The first corner point is the origin. The dim
adjacent points are vectors describing the edges of the parallelepiped with respect to the origin. Additional points are sums of these dim vectors. The variable n_subdivisions
designates the number of subdivisions in each of the dim
directions. Colorizing is done according to hyper_rectangle().
triangulation.generate_fixed_subdivided_parallelepiped(n_subdivisions, corners, colorize = false)
. void GridGenerator::subdivided_parallelepiped | ( | Triangulation< dim > & | tria, |
const unsigned int(&) | n_subdivisions[dim], | ||
const Point< dim >(&) | corners[dim], | ||
const bool | colorize = false |
||
) |
A subdivided parallelepiped, i.e., the same as above, but where the number of subdivisions in each of the dim
directions may vary. Colorizing is done according to hyper_rectangle().
triangulation.generate_varying_subdivided_parallelepiped(n_subdivisions, corners, colorize = false)
. void GridGenerator::subdivided_parallelepiped | ( | Triangulation< dim, spacedim > & | tria, |
const Point< spacedim > & | origin, | ||
const std::array< Tensor< 1, spacedim >, dim > & | edges, | ||
const std::vector< unsigned int > & | subdivisions = {} , |
||
const bool | colorize = false |
||
) |
A subdivided parallelepiped.
tria | The triangulation to be created. It needs to be empty upon calling this function. |
origin | First corner of the parallelepiped. |
edges | An array of dim tensors describing the length and direction of the edges from origin . |
subdivisions | Number of subdivisions in each of the dim directions. Each entry must be positive. An empty vector is equivalent to one subdivision in each direction. |
colorize | Assign different boundary ids if set to true (see the glossary entry on colorization). |
dim
and spacedim
.void GridGenerator::enclosed_hyper_cube | ( | Triangulation< dim > & | tria, |
const double | left = 0. , |
||
const double | right = 1. , |
||
const double | thickness = 1. , |
||
const bool | colorize = false |
||
) |
Hypercube with a layer of hypercubes around it. Parameters left
and right
give the lower and upper bound of the inner hypercube in all coordinate directions. thickness
marks the size of the layer cells.
If the flag colorize
is set, the outer cells get material ids according to the following scheme: extending over the inner cube in (+/-) x-direction 1/2, y-direction 4/8, z-direction 16/32. A bitwise OR operation is used to get these values at the corners and edges (3d), (see also the glossary entry on colorization).
Presently only available in 2d and 3d.
triangulation.generate_hyper_cube(left = 0, right = 1., thickness = 1., colorize = false)
. void GridGenerator::hyper_ball | ( | Triangulation< dim > & | tria, |
const Point< dim > & | center = Point< dim >() , |
||
const double | radius = 1. , |
||
const bool | attach_spherical_manifold_on_boundary_cells = false |
||
) |
Initialize the given triangulation with several coarse mesh cells that cover a hyperball, i.e. a circle in 2d or a ball in 3d, around center
with given radius
. The function is used in step-6.
In order to avoid degenerate cells at the boundaries, the circle is triangulated by five cells, whereas in 3d the ball is subdivided by seven cells. Specifically, these cells are one cell in the center plus one "cap" cell on each of the faces of this center cell. This ensures that under repeated refinement, none of the cells at the outer boundary will degenerate to have an interior angle approaching 180 degrees, as opposed to the case where one might start with just one square (or cube) to approximate the domain. The diameter of the center cell is chosen so that the aspect ratio of the boundary cells after one refinement is optimized.
This function is declared to exist for triangulations of all space dimensions, but throws an error if called in 1d.
By default, the manifold_id is set to 0 on the boundary faces, 1 on the boundary cells, and numbers::flat_manifold_id on the central cell and on internal faces.
A SphericalManifold is attached by default to the boundary faces for correct placement of boundary vertices upon refinement and to be able to use higher order mappings. However, it turns out that this strategy may not be the optimal one to create a good a mesh for a hyperball. The "Possibilities for extensions" section of step-6 has an extensive discussion of how one would construct better meshes and what one needs to do for it. Setting the argument attach_spherical_manifold_on_boundary_cells
to true attaches a SphericalManifold manifold also to the cells adjacent to the boundary, and not only to the boundary faces.
triangulation.generate_hyper_ball(point, radius = 1.)
. void GridGenerator::hyper_ball_balanced | ( | Triangulation< dim > & | tria, |
const Point< dim > & | center = Point< dim >() , |
||
const double | radius = 1. |
||
) |
This is an alternative to hyper_ball with 12 cells in 2d and 32 cells in 3d, which provides a better balance between the size of the cells around the outer curved boundaries and the cell in the interior. The mesh is based on the cells used by GridGenerator::quarter_hyper_ball() with appropriate copies and rotations to fill the whole ball.
The following pictures show the resulting mesh in 2d (left) and 3d:
|
|
By default, the manifold_id is set to 0 on the boundary faces, 1 on the boundary cells, and numbers::flat_manifold_id on the central cell and on internal faces.
triangulation.generate_hyper_ball_balanced(point = Point(), radius = 1.)
. void GridGenerator::non_standard_orientation_mesh | ( | Triangulation< 2 > & | tria, |
const unsigned int | n_rotate_middle_square | ||
) |
Generate a 2d mesh consisting of five squares arranged in a plus-shape. Depending on the number n_rotate_middle_square
passed the middle square is rotated by a degree of n_rotate_middle_square
\(\pi/2\). This way one can generate a mesh in which the middle square contains edges that have the opposite tangential and/or opposite normal orientation compared to the neighboring edges of the other squares.
This mesh is not overly useful from a practical point of view. For debugging purposes it can be used to check for orientation issues for vector- or tensor-valued finite elements.
[out] | tria | The input triangulation. |
[in] | n_rotate_middle_square | number of rotations in [0,4) of right square by \(\pi/2\). |
void GridGenerator::non_standard_orientation_mesh | ( | Triangulation< 3 > & | tria, |
const bool | face_orientation, | ||
const bool | face_flip, | ||
const bool | face_rotation, | ||
const bool | manipulate_left_cube | ||
) |
Generate a 3d mesh consisting of the unit cube joined with a copy shifted by \(s = (1,0,0)\). Depending on the flags passed either the right or the left cube (when looking at the positively oriented (x,z)-plane) contains a face that is either not in standard orientation and/or is rotated by either \(\pi/2\), \(\pi\) or \(3/2\pi\).
This mesh is not overly useful from a practical point of view. For debugging purposes it can be used to check for orientation issues for vector- or tensor-valued finite elements.
[out] | tria | The input triangulation. |
[in] | face_orientation | true if the face is the not in standard orientation. |
[in] | face_flip | true if the face is rotated by +180 degrees |
[in] | face_rotation | true if the face is rotated (additionally) by +90 degrees |
[in] | manipulate_left_cube | true if the left cube is to be re-ordered. If false , it is the right cube. |
void GridGenerator::hyper_sphere | ( | Triangulation< spacedim - 1, spacedim > & | tria, |
const Point< spacedim > & | center = Point< spacedim >() , |
||
const double | radius = 1. |
||
) |
Creates a hyper sphere, i.e., a surface of a ball in spacedim
dimensions. This function only exists for dim+1=spacedim in 2 and 3 space dimensions. (To create a mesh of a ball, use GridGenerator::hyper_ball().)
By default, all manifold ids of the triangulation are set to zero, and a SphericalManifold is attached to the grid.
The following pictures are generated with:
See the documentation topic on manifolds for more details.
triangulation.generate_hyper_sphere(center, radius = 1.)
. void GridGenerator::quarter_hyper_ball | ( | Triangulation< dim > & | tria, |
const Point< dim > & | center = Point< dim >() , |
||
const double | radius = 1. |
||
) |
This function produces a hyper-ball intersected with the positive orthant relative to center
, which contains three elements in 2d and four in 3d. The interior points of the mesh are chosen to balance the minimal singular value of the Jacobian of the mapping from reference to real coordinates among the cells around the interior point, which corresponds to a high mesh quality.
The boundary indicators for the final triangulation are 0 for the curved boundary and 1 for the cut plane. The manifold id for the curved boundary is set to zero, and a SphericalManifold is attached to it.
The resulting grid in 2d and 3d looks as follows:
|
|
triangulation.generate_quarter_hyper_ball(center, radius = 1.)
. void GridGenerator::half_hyper_ball | ( | Triangulation< dim > & | tria, |
const Point< dim > & | center = Point< dim >() , |
||
const double | radius = 1. |
||
) |
This function produces a half hyper-ball around center
, which contains four elements in 2d and 6 in 3d. The cut plane is perpendicular to the x-axis.
The boundary indicators for the final triangulation are 0 for the curved boundary and 1 for the cut plane. The manifold id for the curved boundary is set to zero, and a SphericalManifold is attached to it.
triangulation.generate_half_hyper_ball(center, radius = 1.)
. void GridGenerator::cylinder | ( | Triangulation< dim > & | tria, |
const double | radius = 1. , |
||
const double | half_length = 1. |
||
) |
Create a dim
dimensional cylinder where the \(x\)-axis serves as the axis of the cylinder. For the purposes of this function, a cylinder is defined as a (dim
- 1) dimensional disk of given radius
, extruded along the axis of the cylinder (which is the first coordinate direction). Consequently, in three dimensions, the cylinder extends from x=-half_length
to x=+half_length
and its projection into the yz-plane
is a circle of radius radius
. In two dimensions, the cylinder is a rectangle from x=-half_length
to x=+half_length
and from y=-radius
to y=radius
.
The boundaries are colored according to the following scheme: 0 for the hull of the cylinder, 1 for the left hand face and 2 for the right hand face (see the glossary entry on colorization).
The manifold id for the hull of the cylinder is set to zero, and a CylindricalManifold is attached to it.
triangulation.generate_cylinder(radius = 1., half_length = 1.)
. void GridGenerator::subdivided_cylinder | ( | Triangulation< dim > & | tria, |
const unsigned int | x_subdivisions, | ||
const double | radius = 1. , |
||
const double | half_length = 1. |
||
) |
Create a dim
dimensional cylinder where the \(x\)-axis serves as the axis of the cylinder. For the purposes of this function, a cylinder is defined as a (dim
- 1) dimensional disk of given radius
, extruded along the axis of the cylinder (which is the first coordinate direction). Consequently, in three dimensions, the cylinder extends from x=-half_length
to x=+half_length
and its projection into the yz-plane
is a circle of radius radius
. In two dimensions, the cylinder is a rectangle from x=-half_length
to x=+half_length
and from y=-radius
to y=radius
. This function is only implemented for dim==3.
The boundaries are colored according to the following scheme: 0 for the hull of the cylinder, 1 for the left hand face and 2 for the right hand face (see the glossary entry on colorization).
The manifold id for the hull of the cylinder is set to zero, and a CylindricalManifold is attached to it.
tria | The triangulation to be created. It needs to be empty upon calling this function. |
x_subdivisions | A positive integer denoting the number of cells to generate in the x direction. The default cylinder has x_repetitions=2. |
radius | The radius of the circle in the yz-plane used to extrude the cylinder. |
half_length | The half-length of the cylinder in the x direction. |
triangulation.generate_subdivided_cylinder(x_subdivisions, radius = 1., half_length = 1.)
. void GridGenerator::truncated_cone | ( | Triangulation< dim > & | tria, |
const double | radius_0 = 1.0 , |
||
const double | radius_1 = 0.5 , |
||
const double | half_length = 1.0 |
||
) |
Create a cut cone around the x-axis. The cone extends from x=-half_length
to x=half_length
and its projection into the yz-plane
is a circle of radius radius_0
at x=-half_length
and a circle of radius radius_1
at x=+half_length
. In between the radius is linearly decreasing.
In two dimensions, the cone is a trapezoid from x=-half_length
to x=+half_length
and from y=-radius_0
to y=radius_0
at x=-half_length
and from y=-radius_1
to y=radius_1
at x=+half_length
. In between the range of y
is linearly decreasing.
The boundaries are colored according to the following scheme: 0 for the hull of the cone, 1 for the left hand face, and 2 for the right hand face (see the glossary entry on colorization). Both the boundary indicators and the manifold indicators are set.
In three dimensions, the manifold id of the hull is set to zero, and a CylindricalManifold is attached to it.
Here are the grids in 2d and 3d after two mesh refinements:
triangulation.generate_truncated_cone(radius_0 = 1., radius_2 = 0.5, half_length = 1.)
. void GridGenerator::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 |
||
) |
Initialize the given triangulation with a pipe junction, which is the intersection of three truncated cones.
The geometry has four characteristic cross sections, located at the three openings and the bifurcation. They need to be specified via the function's arguments: each cross section is described by a characteristic point and a radius. The cross sections at the openings are circles and are described by their center point and radius. The bifurcation point describes where the symmetry axes of all cones meet.
Each truncated cone is transformed so that the three merge seamlessly into each other. The bifurcation radius describes the radius that each original, untransformed, truncated cone would have at the bifurcation. This radius is necessary for the construction of the geometry and can, in general, no longer be found in the final result.
Each cone will be assigned a distinct material ID that matches the index of their opening in the argument openings
. For example, the cone which connects to opening with index 0 in openings
will have material ID 0.
Similarly, boundary IDs are assigned to the cross-sections of each opening to match their index. All other boundary faces will be assigned boundary ID 3.
Manifold IDs will be set on the mantles of each truncated cone in the same way. Each cone will have a special manifold object assigned, which is based on the CylindricalManifold class. Further, all cells adjacent to the mantle are given the manifold ID 3. If desired, you can assign an (expensive) TransfiniteInterpolationManifold object to that particular layer of cells with the following code snippet.
dim = 3
and spacedim = 3
.tria | An empty triangulation which will hold the pipe junction geometry. |
openings | Center point and radius of each of the three openings. The container has to be of size three. |
bifurcation | Center point of the bifurcation and hypothetical radius of each truncated cone at the bifurcation. |
aspect_ratio | Aspect ratio of cells, specified as radial over z-extension. Default ratio is \(\Delta r/\Delta z = 1/2\). |
Common configurations of tee fittings (that is, "T" fittings, mimicking the geometry of the letter "T") can be generated with the following sets of parameters:
Corner piece | ||
---|---|---|
Point | Radius | |
Openings | \((2,0,0)\) | \(1\) |
\((0,2,0)\) | \(1\) | |
\((0,0,2)\) | \(1\) | |
Bifurcation | \((0,0,0)\) | \(1\) |
T-pipe | ||
---|---|---|
Point | Radius | |
Openings | \((-2,0,0)\) | \(1\) |
\((0,2,0)\) | \(1\) | |
\((2,0,0)\) | \(1\) | |
Bifurcation | \((0,0,0)\) | \(1\) |
Y-pipe | ||
---|---|---|
Point | Radius | |
Openings | \((-2,0,0)\) | \(1\) |
\((1,\sqrt{3},0)\) | \(1\) | |
\((1,-\sqrt{3},0)\) | \(1\) | |
Bifurcation | \((0,0,0)\) | \(1\) |
Definition at line 266 of file grid_generator_pipe_junction.cc.
void GridGenerator::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.
Each of the square mesh cells is Cartesian and has size one in each coordinate direction. The center of cell number zero is the origin.
tria | A Triangulation object which has to be empty. |
sizes | A vector of integers of dimension GeometryInfo<dim>::faces_per_cell with the following meaning: the legs of the cross are stacked on the faces of the center cell, in the usual order of deal.II cells, namely first \(-x\), then \(x\), then \(-y\) and so on. The corresponding entries in sizes name the number of cells stacked on this face. All numbers may be zero, thus L- and T-shaped domains are specializations of this domain. |
colorize_cells | If colorization is enabled, then the material id of a cells corresponds to the leg it is in. The id of the center cell is zero, and then the legs are numbered starting at one (see the glossary entry on colorization). |
Examples in two and three dimensions are:
void GridGenerator::hyper_L | ( | Triangulation< dim > & | tria, |
const double | left = -1. , |
||
const double | right = 1. , |
||
const bool | colorize = false |
||
) |
Initialize the given triangulation with a hyper-L (in 2d or 3d) consisting of exactly 2^dim-1
cells. It produces the hypercube with the interval [left,right] without the hypercube made out of the interval [(left+right)/2,right] for each coordinate. Because the domain is about the simplest one with a reentrant (i.e., non-convex) corner, solutions of many partial differential equations have singularities at this corner. That is, at the corner, the gradient or a higher derivative (depending on the boundary conditions chosen) does not remain bounded. As a consequence, this domain is often used to test convergence of schemes when the solution lacks regularity.
If the colorize
flag is true
, the boundary_ids
of the surfaces are assigned such that the left boundary is 0 and the others are assigned counterclockwise in ascending order (see the glossary entry on colorization). The colorize
option only works in two dimensions.
This function will create the classical L-shape in 2d and it will look like the following in 3d:
This function exists for triangulations of all space dimensions, but throws an error if called in 1d.
void GridGenerator::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 | ||
) |
Initialize the given triangulation in 2d or 3d with a generalized subdivided hyper-L.
This function produces a subdivided hyper rectangle with dimensions given by bottom_left
and top_right
, with the given number of subdivisions in each direction given in the vector repetitions
, and with a number of cells removed, given in the vector n_cells_to_remove
. Note that n_cells_to_remove
contains integers, meaning that its entries can be both positive and negative. A positive number denotes cutting away cells in the 'positive' orientation, for example left to right in the x-direction, bottom to top in the y-direction, and front to back in the z-direction. A negative number denotes cutting away cells in the reverse direction, so right to left, top to bottom, and back to front.
A demonstration of this grid can be found in step-75.
This function may be used to generate a mesh for a backward facing step, a useful domain for benchmark problems in fluid dynamics. The first image is a backward facing step in 3d, generated by removing all cells in the z-direction, and 2 cells in the positive x- and y-directions:
And in 2d, we can cut away 1 cell in the negative x-direction, and 2 cells in the negative y-direction:
void GridGenerator::hyper_cube_slit | ( | Triangulation< dim > & | tria, |
const double | left = 0. , |
||
const double | right = 1. , |
||
const bool | colorize = false |
||
) |
Initialize the given Triangulation with a hypercube with a slit. In each coordinate direction, the hypercube extends from left
to right
.
In 2d, the split goes in vertical direction from x=(left+right)/2, y=left
to the center of the square at x=y=(left+right)/2
.
In 3d, the 2d domain is just extended in the z-direction, such that a plane cuts the lower half of a rectangle in two. This function is declared to exist for triangulations of all space dimensions, but throws an error if called in 1d.
If colorize
is set to true
, the faces forming the slit are marked with boundary id 1 and 2, respectively (see the glossary entry on colorization).
void GridGenerator::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 |
||
) |
Produce a hyper-shell, the region between two spheres around center
, with given inner_radius
and outer_radius
. The number n_cells
indicates the number of cells of the resulting triangulation, i.e., how many cells form the ring (in 2d) or the shell (in 3d).
If the flag colorize
is true
, then the outer boundary will have the indicator 1 while the inner boundary has id zero. In 3d, this applies to both the faces and the edges of these boundaries. If the flag is false
, both have indicator zero (see the glossary entry on colorization).
All manifold ids are set to zero, and a SphericalManifold is attached to every cell and face of the triangulation.
In 2d, the number n_cells
of elements for this initial triangulation can be chosen arbitrarily. If the number of initial cells is zero (as is the default), then it is computed adaptively such that the resulting elements have the least aspect ratio.
In 3d, only certain numbers are allowed:
The versions with 24, 48, and \(2^m 192\) cells are useful if the shell is thin and the radial lengths should be made more similar to the circumferential lengths.
The 3d grids with 12 and 96 cells are plotted below:
triangulation.generate_hyper_shell(center, inner_radius, outer_radius, n_cells = 0, colorize = false)
. void GridGenerator::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 | ||
) |
Produce an eccentric hyper-shell, the region between two spheres centered on two distinct center points. One has to specify the inner_center
and outer_center
, with given inner_radius
and outer_radius
. The number n_cells
indicates the number of cells of the resulting triangulation, i.e., how many cells form the ring (in 2d) or the shell (in 3d).
By default, the outer boundary has the indicator 1 while the inner boundary has id 0. In 3d, this applies to both the faces and the edges of these boundaries.
A SphericalManifold is attached to the outer boundary with an id of 1 while another SphericalManifold is attached to the inner boundary with an id of 0. A TransfiniteInterpolationManifold is attached to all other cells and faces of the triangulation with an id of 2.
Here, the number n_cells
of elements has the same meaning as in GridGenerator::hyper_shell.
The grids with a 30% offset of the inner shell in the x direction, 12 initial cells and 3 levels of global refinement are plotted below:
void GridGenerator::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 |
||
) |
Produce a half hyper-shell, i.e. the space between two circles in two space dimensions and the region between two spheres in 3d, with given inner and outer radius and a given number of elements for this initial triangulation. However, opposed to the previous function, it does not produce a whole shell, but only one half of it, namely that part for which the first component is restricted to non-negative values. The purpose of this function is to enable computations for solutions which have rotational symmetry, in which case the half shell in 2d represents a shell in 3d.
If the number of initial cells n_cells
is zero in 2d (as is the default), then it is computed adaptively such that the resulting elements have the least aspect ratio. The argument is ignored in 3d, where the coarse mesh always has 5 cells.
If colorize is set to true
, the inner, outer, and the part of the boundary where \(x=0\), get indicator 0, 1, and 2, respectively. Additionally, in 2d, the boundary indicator 3 is given to the vertical edge below the x-axis. Otherwise, if colorize is set to false
all indicators are set to 0 (see the glossary entry on colorization).
All manifold ids are set to zero, and a SphericalManifold is attached to the triangulation.
void GridGenerator::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 |
||
) |
Produce a domain that is the intersection between a hyper-shell with given inner and outer radius, i.e. the space between two circles in two space dimensions and the region between two spheres in 3d, and the positive quadrant (in 2d) or octant (in 3d). In 2d, this is indeed a quarter of the full annulus, while the function is a misnomer in 3d because there the domain is not a quarter but one eighth of the full shell.
If the number of initial cells is zero (as is the default), then it is computed adaptively such that the resulting elements have the least aspect ratio in 2d.
If colorize
is set to true
, the inner, outer, left, and right boundary get indicator 0, 1, 2, and 3 in 2d, respectively. Otherwise all indicators are set to 0. In 3d indicator 2 is at the face \(x=0\), 3 at \(y=0\), 4 at \(z=0\) (see the glossary entry on colorization).
All manifold ids are set to zero, and a SphericalManifold is attached to the triangulation.
void GridGenerator::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 |
||
) |
Produce a domain that is the space between two cylinders in 3d, with given length, inner and outer radius and a given number of elements. The cylinder shell is built around the \(z\)-axis with the two faces located at \(z = 0\) and \(z = \) length
.
If n_radial_cells
is zero (as is the default), then it is computed adaptively such that the resulting elements have the least aspect ratio. The same holds for n_axial_cells
.
If colorize
is set to true, a boundary id of 0 is set for the inner cylinder, a boundary id of 1 is set for the outer cylinder, a boundary id of 2 is set for the bottom (z-) boundary and a boundary id of 3 is set for the top (z+) boundary.
All manifold ids are set to zero, and a CylindricalManifold is attached to the triangulation.
In this picture, a cylinder shell of length 2, inner radius 0.5, outer radius 1 is shown. The default argument for n_radial_cells and n_axial_cells are used and a single global refinement is carried out.
void GridGenerator::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 |
||
) |
Produce the volume or surface mesh of a torus. The axis of the torus is the \(y\)-axis while the plane of the torus is the \(x\)- \(z\) plane.
If dim
is 3, the mesh will be the volume of the torus, using a mesh equivalent to the circle in the poloidal coordinates with 5 cells on the cross section. This function attaches a TorusManifold to all boundary faces which are marked with a manifold id of 1, a CylindricalManifold to the interior cells and all their faces which are marked with a manifold id of 2 (representing a flat state within the poloidal coordinates), and a TransfiniteInterpolationManifold to the cells between the TorusManifold on the surface and the ToroidalManifold in the center, with cells marked with manifold id 0.
An example for the case if dim
is 3 with a cut through the domain at \(z=0\), 6 toroidal cells, \(R=2\) and \(r=0.5\) without any global refinement is given here:
In this picture, the light gray shade represents the manifold id 0 of the transfinite interpolation, which is applied to smoothly add new points between the toroidal shape on the domain boundary and the inner rim where a cylindrical description around the y-axis is prescribed. The inner rim with the manifold id 2 is shown in red shade.
If dim
is 2, the mesh will describe the surface of the torus and this function attaches a TorusManifold to all cells and faces (which are marked with a manifold id of 0).
tria | The triangulation to be filled. |
centerline_radius | The radius of the circle which forms the center line of the torus containing the loop of cells. Must be greater than inner_radius . |
inner_radius | The distance between the inner edge of the torus and origin. |
n_cells_toroidal | Optional argument to set the number of cell layers in toroidal direction. The default is 6 cell layers. |
phi | Optional argument to generate an open torus with angle \(0 < \varphi <= 2 \pi\). The default value is \(2 \pi\), in which case a closed torus is generated. If the torus is open, the torus is cut at two planes perpendicular to the torus centerline. The center of these two planes are located at \((x_1, y_1, z_1) = (R, 0, 0)\) and \((x_2, y_2, z_2) = (R \cos(\varphi), 0, R \sin(\varphi))\). |
void GridGenerator::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 |
||
) |
This function produces a square in the xy-plane with a cylindrical hole in the middle. The square and the circle are centered at the origin. In 3d, this geometry is extruded in \(z\) direction to the interval \([0,L]\).
The inner boundary has a manifold id of \(0\) and a boundary id of \(6\). This function attaches a PolarManifold or CylindricalManifold to the interior boundary in 2d and 3d respectively. The other faces have boundary ids of \(0, 1, 2, 3, 4\), or \(5\) given in the standard order of faces in 2d or 3d.
It is implemented in 2d and 3d, and takes the following arguments:
triangulation | The triangulation to be filled. |
inner_radius | Radius of the internal hole. |
outer_radius | Half of the edge length of the square. |
L | Extension in z-direction (only used in 3d). |
repetitions | Number of subdivisions along the z-direction . |
colorize | Whether to assign different boundary indicators to different faces (see the glossary entry on colorization). The colors are given in lexicographic ordering for the flat faces (0 to 3 in 2d, 0 to 5 in 3d) plus the curved hole (4 in 2d, and 6 in 3d). If colorize is set to false, then flat faces get the number 0 and the hole gets number 1. |
triangulation.generate_hyper_cube_with_cylindrical_hole(inner_radius = 0.25, outer_radius = 0.5, L = 0.5, repetitions = 1, colorize = false)
. void GridGenerator::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 |
||
) |
Produce a grid consisting of concentric shells. The primary difference between this function and GridGenerator::hyper_shell() is that this function permits unevenly spaced (in the radial direction) coarse level cells.
The parameters center
, inner_radius
, and outer_radius
behave in the same way as the first three arguments to GridGenerator::hyper_shell. n_shells
gives the total number of shells to use (i.e., the number of cells in the radial direction). The outer radius of the \(k\)th shell is given by
\[ r = r_{\mathrm{inner}} + (r_\mathrm{outer} - r_\mathrm{inner}) \left(1 - \frac{ \tanh(\mathrm{skewness}(1 - k/\mathrm{n\_shells}))} {\tanh(\mathrm{skewness})}\right) \]
where skewness
is a parameter controlling the shell spacing in the radial direction: values of skewness
close to zero correspond to even spacing, while larger values of skewness
(such as \(2\) or \(3\)) correspond to shells biased to the inner radius.
n_cells_per_shell
is the same as in GridGenerator::hyper_shell: in 2d the default choice of zero will result in 8 cells per shell (and 12 in 3d). The only valid values in 3d are 6 (the default), 12, and 96 cells: see the documentation of GridGenerator::hyper_shell for more information.
If colorize
is true
then the outer boundary of the merged shells has a boundary id of \(1\) and the inner boundary has a boundary id of \(0\).
Example: The following code (see, e.g., step-10 for instructions on how to visualize GNUPLOT output)
generates the following output:
void GridGenerator::moebius | ( | Triangulation< 3, 3 > & | tria, |
const unsigned int | n_cells, | ||
const unsigned int | n_rotations, | ||
const double | R, | ||
const double | r | ||
) |
Produce a ring of cells in 3d that is cut open, twisted and glued together again. This results in a kind of moebius-loop.
tria | The triangulation to be worked on. |
n_cells | The number of cells in the loop. Must be greater than 4. |
n_rotations | The number of rotations ( \(\pi/2\) each) to be performed before gluing the loop together. |
R | The radius of the circle, which forms the middle line of the torus containing the loop of cells. Must be greater than r . |
r | The radius of the cylinder bent together as a loop. |
void GridGenerator::generate_from_name_and_arguments | ( | Triangulation< dim, spacedim > & | tria, |
const std::string & | grid_generator_function_name, | ||
const std::string & | grid_generator_function_arguments | ||
) |
Call one of the other GridGenerator functions, parsing the name of the function to call from the string grid_generator_function_name
, and the arguments to the function from the string grid_generator_function_arguments
.
The string that supplies the arguments is passed to the function Patterns::Tools::Convert<TupleType>::to_value(), where TupleType
here is a tuple containing all the arguments of the GridGenerator function, including all optional arguments.
An example usage of this function is given by:
Here, the colon separates the function arguments, and the comma separates the coordinates of a Point<2> argument.
According to the arity of the TupleType
, the arguments of the function may be separated by different separators (see the documentation of Patterns::Tuple for the details of how the conversion is performed). If a wrong format is used, an exception is thrown, and the expected format is output as an error message.
All GridGenerator functions are supported. If you find some that are missing, please open an issue on GitHub.
tria | The triangulation to be worked on |
grid_generator_function_name | The name of the function to call |
grid_generator_function_arguments | The arguments of the function, in the format of a tuple-convertible string |
Definition at line 352 of file grid_generator_from_name.cc.
void GridGenerator::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 |
||
) |
Generate a Triangulation from the zero level set of an implicit function, using the CGAL library.
This function is only implemented for dim
equal to two or three, and requires that deal.II is configured using DEAL_II_WITH_CGAL
. When dim
is equal to three, the implicit_function
is supposed to be negative in the interior of the domain, positive outside, and to be entirely enclosed in a ball of radius outer_ball_radius
centered at the point interior_point
. The triangulation that is generated covers the volume bounded by the zero level set of the implicit function where the implicit_function
is negative.
When dim
is equal to two, the generated surface triangulation is the zero level set of the implicit_function
, oriented such that the surface triangulation has normals pointing towards the region where implicit_function
is positive.
The struct data
can be used to pass additional arguments to the CGAL::Mesh_criteria_3 class (see https://doc.cgal.org/latest/Mesh_3/index.html for more information.)
An example usage of this function is given by
The above snippet of code generates the following grid for dim
equal to two and three respectively
Also see Simplex support.
[out] | tria | The output triangulation |
[in] | implicit_function | The implicit function |
[in] | data | Additional parameters to pass to the CGAL::make_mesh_3 function and to the CGAL::make_surface_mesh functions |
[in] | interior_point | A point in the interior of the domain, for which implicit_function is negative |
[in] | outer_ball_radius | The radius of the ball that will contain the generated Triangulation object |
void GridGenerator::surface_mesh_to_volumetric_mesh | ( | const Triangulation< 2, 3 > & | surface_tria, |
Triangulation< 3 > & | vol_tria, | ||
const CGALWrappers::AdditionalData< 3 > & | data = CGALWrappers::AdditionalData< 3 >{} |
||
) |
Create a deal.II Triangulation<3> out of a deal.II Triangulation<2,3> by filling it with tetrahedra.
The last optional argument data
can be used to pass additional arguments to the CGAL::Mesh_criteria_3 class (see https://doc.cgal.org/latest/Mesh_3/index.html for more information).
[in] | surface_tria | The input deal.II Triangulation<2,3>. |
[out] | vol_tria | The output deal.II Triangulation<3>. |
[in] | data | Additional parameters to pass to the CGAL::make_mesh_3 function. |
void GridGenerator::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 |
||
) |
Given the two triangulations specified as the first two arguments, create the triangulation that contains the cells of both triangulation and store it in the third parameter. Previous content of result
will be deleted. One of the two input triangulations can also be the result
triangulation.
This function is most often used to compose meshes for more complicated geometries if the geometry can be composed of simpler parts for which functions exist to generate coarse meshes. For example, the channel mesh used in step-35 could in principle be created using a mesh created by the GridGenerator::hyper_cube_with_cylindrical_hole function and several rectangles, and merging them using the current function. The rectangles will have to be translated to the right for this, a task that can be done using the GridTools::shift function (other tools to transform individual mesh building blocks are GridTools::transform, GridTools::rotate, and GridTools::scale).
Vertices that are less than duplicated_vertex_tolerance
apart will be merged together. It is usually necessary to set this value to something that depends on the input triangulations in some way. One reasonable choice is to use the minimum distance between all adjacent vertices of the input mesh divided by some constant:
This will merge any vertices that are closer than any pair of vertices on the input meshes. If the tolerance is set to zero, vertices are not merged.
copy_manifold_ids
is set to true
, manifold ids will be copied. If copy_boundary_ids
is set to true
, boundary_ids are copied to all remaining faces at the boundary.result
, nor does it set any manifold ids. In particular, manifolds attached to the two input triangulations will be lost in the result
triangulation.void GridGenerator::merge_triangulations | ( | const std::vector< const Triangulation< dim, spacedim > * > & | triangulations, |
Triangulation< dim, spacedim > & | result, | ||
const double | duplicated_vertex_tolerance = 1.0e-12 , |
||
const bool | copy_manifold_ids = false , |
||
const bool | copy_boundary_ids = false |
||
) |
Same as above but allows to merge more than two triangulations at once. The following gives an example of how to use this function:
result.merge_triangulations(triangulations, duplicated_vertex_tolerance = 1e-12, copy_manifold_ids = false)
. void GridGenerator::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.
input | The triangulation which will be replicated along the coordinate axes. |
extents | A vector with dim entries specifying how many copies of a triangulation should be present along each coordinate axis. |
result | The triangulation to be created. It needs to be empty upon calling this function. |
This function creates a new Triangulation equal to a dim
-dimensional array of copies of input
. Copies of input
are created by translating input
along the coordinate axes. Boundary ids of faces (but not lines in 3d) and all manifold ids are copied but Manifold objects are not since most Manifold objects do not work correctly when a Triangulation has been translated.
To see how this works, consider the following code:
results in
And, similarly, in 3d:
results in
input
based on the BoundingBox of input
. If the boundary faces of input
are not aligned with the coordinate axes then the copies might not share common faces; i.e., this function is intended for simple geometries with boundary faces aligned along the coordinate axes.result.replicate_triangulation(input, extents)
. void GridGenerator::create_union_triangulation | ( | const Triangulation< dim, spacedim > & | triangulation_1, |
const Triangulation< dim, spacedim > & | triangulation_2, | ||
Triangulation< dim, spacedim > & | result | ||
) |
Given the two triangulations specified as the first two arguments, create the triangulation that contains the finest cells of both triangulation and store it in the third parameter. Previous content of result
will be deleted.
triangulation_1
and triangulation_2
have the same manifold descriptions. The output Triangulation has
the same manifold ids as these two triangulations.void GridGenerator::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 | ||
) |
This function creates a triangulation that consists of the same cells as are present in the first argument, except those cells that are listed in the second argument. The purpose of the function is to generate geometries subtractively from the geometry described by an existing triangulation. A prototypical case is a 2d domain with rectangular holes. This can be achieved by first meshing the entire domain and then using this function to get rid of the cells that are located at the holes. A demonstration of this particular use case is part of step-27. Likewise, you could create the mesh that GridGenerator::hyper_L() produces by starting with a GridGenerator::hyper_cube(), refining it once, and then calling the current function with a single cell in the second argument.
[in] | input_triangulation | The original triangulation that serves as the template from which the new one is to be created. |
[in] | cells_to_remove | A list of cells of the triangulation provided as first argument that should be removed (i.e., that should not show up in the result. |
[out] | result | The resulting triangulation that consists of the same cells as are in input_triangulation , with the exception of the cells listed in cells_to_remove . |
result
, nor does it set any manifold ids.void GridGenerator::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 = {} |
||
) |
Extrude the Triangulation input
in the \(z\) direction from \(z = 0\) to \(z =
\text{height}\) and store it in result
. This is done by replicating the input triangulation n_slices
times in \(z\) direction, and then forming (n_slices-1)
layers of cells out of these replicates.
The boundary indicators of the faces of input
will be assigned to the corresponding side walls in \(z\) direction. The bottom and top get the next two free boundary indicators: i.e., if input
has boundary ids of \(0\), \(1\), and \(42\), then the \(z = 0\) boundary id of result
will be \(43\) and the \(z = \text{height}\) boundary id will be \(44\).
This function does not, by default, copy manifold ids. The reason for this is that there is no way to set the manifold ids on the lines of the resulting Triangulation without more information: for example, if two faces of input
with different manifold ids meet at a shared vertex then there is no a priori reason to pick one manifold id or another for the lines created in result
that are parallel to the \(z\)-axis and pass through that point. If copy_manifold_ids
is true
then this function sets line manifold ids by picking the one that appears first in manifold_priorities
. For example: if manifold_priorities
is {0, 42, numbers::flat_manifold_id}
and the line under consideration is adjacent to faces with manifold ids of 0
and 42
, then that line will have a manifold id of 0
. The correct ordering is almost always
In particular, since TransfiniteInterpolationManifold interpolates between surrounding manifolds, its manifold id should usually not be set on lines or faces that are adjacent to cells with different manifold ids. The default value for manifold_priorities
follows this ranking (where each category is sorted in ascending order):
Note that numbers::flat_manifold_id (should it be a manifold id of input
) will always be the last entry in the first category.
[in] | input | A two-dimensional input triangulation. |
[in] | n_slices | The number of times the input triangulation will be replicated in \(z\) direction. These slices will then be connected into (n_slices-1) layers of three-dimensional cells. Clearly, n_slices must be at least two. |
[in] | height | The distance in \(z\) direction between the individual slices. |
[out] | result | The resulting three-dimensional triangulation. |
[in] | copy_manifold_ids | See the description above. |
[in] | manifold_priorities | See the description above. |
input
must be a coarse mesh, i.e., it cannot have any refined cells.input
and output
have different spatial dimensions, no manifold objects are copied by this function regardless of the value of copy_manifold_ids
.input.extrude_triangulation(n_slices, height, result)
. void GridGenerator::extrude_triangulation | ( | const Triangulation< 2, 2 > & | input, |
const unsigned int | n_slices, | ||
const double | height, | ||
Triangulation< 2, 2 > & | result, | ||
const bool | copy_manifold_ids = false , |
||
const std::vector< types::manifold_id > & | manifold_priorities = {} |
||
) |
Overload of extrude_triangulation() to allow dimension independent code to compile. This function throws an error when called, as extrude_triangulation() is only implemented to extrude a dim=2 to a dim=3 Triangulation.
void GridGenerator::extrude_triangulation | ( | const Triangulation< 2, 2 > & | input, |
const std::vector< double > & | slice_coordinates, | ||
Triangulation< 3, 3 > & | result, | ||
const bool | copy_manifold_ids = false , |
||
const std::vector< types::manifold_id > & | manifold_priorities = {} |
||
) |
Overload of the previous function. Take a 2d Triangulation that is being extruded. Differing from the previous function taking height and number of slices for uniform extrusion, this function takes z-axis values slice_coordinates
where the slicing will happen. The boundary indicators of the faces of input
are going to be assigned to the corresponding side walls in z direction. The bottom and top get the next two free boundary indicators.
input
must be a coarse mesh, i.e., it cannot have any refined cells.input
and output
have different spatial dimensions no manifold objects are copied (nor are any manifold ids set) by this function. void GridGenerator::extrude_triangulation | ( | const Triangulation< 2, 2 > & | input, |
const std::vector< double > & | slice_coordinates, | ||
Triangulation< 2, 2 > & | result, | ||
const bool | copy_manifold_ids = false , |
||
const std::vector< types::manifold_id > & | manifold_priorities = {} |
||
) |
Overload of extrude_triangulation() to allow dimension independent code to compile. This function throws an error when called, as extrude_triangulation() is only implemented to extrude a dim=2 to a dim=3 Triangulation.
void GridGenerator::flatten_triangulation | ( | const Triangulation< dim, spacedim1 > & | in_tria, |
Triangulation< dim, spacedim2 > & | out_tria | ||
) |
Given an input triangulation in_tria
, this function makes a new flat triangulation out_tria
which contains a single level with all active cells of the input triangulation. If spacedim1
and spacedim2
are different, only the first few components of the vertex coordinates are copied over. This is useful to create a Triangulation<2,3> out of a Triangulation<2,2>, or to project a Triangulation<2,3> into a Triangulation<2,2>, by neglecting the \(z\) components of the vertices.
No internal checks are performed on the vertices, which are assumed to make sense topologically in the target spacedim2
dimensional space. If this is not the case, you will encounter problems when using the triangulation later on.
All information about cell manifold indicators and material indicators are copied from one triangulation to the other. The same is true for the manifold indicators and, if an object is at the boundary, the boundary indicators of faces and edges of the triangulation.
This function will fail if the input Triangulation is of type parallel::distributed::Triangulation, as well as when the input Triangulation contains hanging nodes. In other words, this function only works for globally refined triangulations.
[in] | in_tria | The base input for a new flat triangulation. |
[out] | out_tria | The desired flattened triangulation constructed from the in_tria. |
input
and output
have different spatial dimensions, no manifold objects are copied by this function: you must attach new manifold objects to out_tria
.in_tria.flatten_triangulation(out_tria)
. void GridGenerator::convert_hypercube_to_simplex_mesh | ( | const Triangulation< dim, spacedim > & | in_tria, |
Triangulation< dim, spacedim > & | out_tria | ||
) |
Convert a triangulation consisting only of hypercube cells (quadrilaterals, hexahedra) to a triangulation only consisting of simplices (triangles, tetrahedra).
As an example, the following image shows how a set of four hexahedra meshing one eighth of a sphere are subdivided into tetrahedra, and how the curved surface is taken into account. Colors indicate how boundary indicators are inherited:
In general, each quadrilateral in 2d is subdivided into eight triangles, and each hexahedron in 3d into 24 tetrahedra, as shown here (top left for the 2d case, the rest shows vertex numbers and subdivisions for a single 3d hexahedron):
Material ID and boundary IDs are inherited upon conversion.
[in] | in_tria | The triangulation containing quadrilateral or hexahedral elements. |
[out] | out_tria | The converted triangulation containing triangular or tetrahedral elements. |
in_tria
to out_tria
, e.g., with the following code: Also see Simplex support.
in_tria.convert_hypercube_to_simplex_mesh(out_tria)
. void GridGenerator::alfeld_split_of_simplex_mesh | ( | const Triangulation< dim, spacedim > & | in_tria, |
Triangulation< dim, spacedim > & | out_tria | ||
) |
Perform an Alfeld split (also called barycentric refinement) of a simplex mesh.
Each simplex cell in the input mesh (given in in_tria
) is refined into three (for dim
= 2) or four (for dim
= 3) simplices connecting to the barycenter, which is the only new vertex added for each input cell. In the process, the simplex mesh is flattened (no hierarchy is kept).
dim
= 2 and hanging nodes are not supported.The meshes produced by this function can be used for Scott-Vogelius elements for the Stokes equation: The \(P_k - DGP_{k-1}\) element is point-wise divergence free on barycentric refined meshes for \(k\geq 2\) for dim
= 2 and \(k\geq 3\) for dim
= 3, see [84].
Also see Simplex support.
void GridGenerator::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 |
||
) |
Create a coordinate-parallel brick from the two diagonally opposite corner points p1
and p2
and subdivide each cell into simplices.
The number of vertices in coordinate direction i
is given by repetitions[i]+1
.
dim
2) or 5 tetrahedra (for dim
3), respectively.dim==spacedim
.Also see Simplex support.
void GridGenerator::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 |
||
) |
Initialize the given triangulation with a hypercube (square in 2d and cube in 3d) consisting of repetitions
cells in each direction with each cell divided into simplices.
The hypercube volume is the tensor product interval \([left,right]^{\text{dim}}\) in the present number of dimensions, where the limits are given as arguments. They default to zero and unity, then producing the unit hypercube.
dim
2) or 5 tetrahedra (for dim
3), respectively.Also see Simplex support.
return_type GridGenerator::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 >() |
||
) |
This function implements a boundary subgrid extraction. Given a <dim,spacedim>-Triangulation (the "volume mesh") the function extracts a subset of its boundary (the "surface mesh"). The boundary to be extracted is specified by a list of boundary_ids. If none is specified the whole boundary will be extracted. The function is used in step-38.
The function also builds a mapping linking the cells on the surface mesh to the corresponding faces on the volume one. This mapping is the return value of the function.
MeshType | A type that satisfies the requirements of the MeshType concept. The map that is returned will be between cell iterators pointing into the container describing the surface mesh and face iterators of the volume mesh container. If MeshType is DoFHandler, then the function will re-build the triangulation underlying the second argument and return a map between appropriate iterators into the MeshType arguments. However, the function will not actually distribute degrees of freedom on this newly created surface mesh. |
dim | The dimension of the cells of the volume mesh. For example, if dim==2, then the cells are quadrilaterals that either live in the plane, or form a surface in a higher-dimensional space. The dimension of the cells of the surface mesh is consequently dim-1. |
spacedim | The dimension of the space in which both the volume and the surface mesh live. |
[in] | volume_mesh | A container of cells that define the volume mesh. |
[out] | surface_mesh | A container whose associated triangulation will be built to consist of the cells that correspond to the (selected portion of) the boundary of the volume mesh. |
[in] | boundary_ids | A list of boundary indicators denoting that subset of faces of volume cells for which this function should extract the surface mesh. If left at its default, i.e., if the set is empty, then the function operates on all boundary faces. |
To prevent printing a very long type in the doxygen documentation the actual return type of this function is
when MeshType
is DoFHandler and
when MeshType
is Triangulation and and not the shorter stub provided here.
volume_mesh
and surface_mesh
have different spatial dimensions no manifold objects are copied by this function: you must attach new manifold objects to surface_mesh
.