Classes | |
class | ExcInvalidRadii |
class | ExcInvalidRepetitions |
class | ExcInvalidRepetitionsDimension |
Static Public Member Functions | |
template<int dim, int spacedim> | |
static void | hyper_cube (Triangulation< dim, spacedim > &tria, const double left=0., const double right=1.) |
template<int dim> | |
static void | subdivided_hyper_cube (Triangulation< dim > &tria, const unsigned int repetitions, const double left=0., const double right=1.) |
template<int dim, int spacedim> | |
static void | hyper_rectangle (Triangulation< dim, spacedim > &tria, const Point< spacedim > &p1, const Point< spacedim > &p2, const bool colorize=false) |
template<int dim> | |
static void | subdivided_hyper_rectangle (Triangulation< dim > &tria, const std::vector< unsigned int > &repetitions, const Point< dim > &p1, const Point< dim > &p2, const bool colorize=false) |
template<int dim> | |
static 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) |
template<int dim> | |
static void | subdivided_hyper_rectangle (Triangulation< dim > &tria, const std::vector< std::vector< double > > &spacing, const Point< dim > &p, const Table< dim, unsigned char > &material_id, const bool colorize=false) |
template<int dim> | |
static void | parallelogram (Triangulation< dim > &tria, const Tensor< 2, dim > &corners, const bool colorize=false) |
template<int dim> | |
static 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> | |
static void | hyper_ball (Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.) |
template<int dim> | |
static void | half_hyper_ball (Triangulation< dim > &tria, const Point< dim > ¢er=Point< dim >(), const double radius=1.) |
template<int dim> | |
static void | cylinder (Triangulation< dim > &tria, const double radius=1., const double half_length=1.) |
template<int dim> | |
static void | hyper_L (Triangulation< dim > &tria, const double left=-1., const double right=1.) |
template<int dim> | |
static void | hyper_cube_slit (Triangulation< dim > &tria, const double left=0., const double right=1., const bool colorize=false) |
template<int dim> | |
static void | hyper_shell (Triangulation< dim > &tria, const Point< dim > ¢er, const double inner_radius, const double outer_radius, const unsigned int n_cells=0, bool colorize=false) |
template<int dim> | |
static 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) |
template<int dim> | |
static 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) |
template<int dim> | |
static 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 repetition=1, const bool colorize=false) |
static 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> | |
static void | laplace_transformation (Triangulation< dim > &tria, const std::map< unsigned int, Point< dim > > &new_points) |
Static Private Member Functions | |
template<int dim, int spacedim> | |
static void | colorize_hyper_rectangle (Triangulation< dim, spacedim > &tria) |
template<int dim> | |
static void | colorize_subdivided_hyper_rectangle (Triangulation< dim > &tria, const Point< dim > &p1, const Point< dim > &p2, const double epsilon) |
template<int dim> | |
static void | colorize_hyper_shell (Triangulation< dim > &tria, const Point< dim > ¢er, const double inner_radius, const double outer_radius) |
static void | laplace_solve (const SparseMatrix< double > &S, const std::map< unsigned int, double > &m, Vector< double > &u) |
This class provides a collection of functions for generating basic triangulations. Below, we try to provide some pictures in order to illustrate at least the more complex ones.
Some of these functions receive a flag colorize
. If this is set, parts of the boundary receive different boundary numbers, allowing them to be distinguished by application programs. See the documentation of the functions for details.
Additionally this class provides a function (laplace_transformation
) that smoothly transforms a grid according to given new boundary points. This can be used to transform (simple-shaped) grids to a more complicated ones, like a shell onto a grid of an airfoil, for example.
No meshes for the codimension one case are provided at the moment.
static void GridGenerator::hyper_cube | ( | Triangulation< dim, spacedim > & | tria, | |
const double | left = 0. , |
|||
const double | right = 1. | |||
) | [inline, static] |
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]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.
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.
static void GridGenerator::subdivided_hyper_cube | ( | Triangulation< dim > & | tria, | |
const unsigned int | repetitions, | |||
const double | left = 0. , |
|||
const double | right = 1. | |||
) | [inline, static] |
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 spacedim=dim+1 the same mesh as in the case spacedim=dim is created, but the vertices have an additional coordinate =0. So, if dim=1 one obtains line along the x axis in the xy plane, and if dim=3 one obtains a square in lying in the xy plane in 3d space.
static void GridGenerator::hyper_rectangle | ( | Triangulation< dim, spacedim > & | tria, | |
const Point< spacedim > & | p1, | |||
const Point< spacedim > & | p2, | |||
const bool | colorize = false | |||
) | [inline, static] |
Create a coordinate-parallel brick from the two diagonally opposite corner points p1
and p2
.
If the colorize
flag is set, the boundary_indicators
of the surfaces 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. 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. For instance, the center point (1,-1,1) yields a material id 5.
static void GridGenerator::subdivided_hyper_rectangle | ( | Triangulation< dim > & | tria, | |
const std::vector< unsigned int > & | repetitions, | |||
const Point< dim > & | p1, | |||
const Point< dim > & | p2, | |||
const bool | colorize = false | |||
) | [inline, static] |
Create a coordinate-parallel parallelepiped from the two diagonally opposite corner points p1
and p2
. In dimension i
, repetitions[i]
cells are generated.
To get cells with an aspect ratio different from that of the domain, use different numbers of subdivisions in different coordinate directions. The minimum number of subdivisions in each direction is 1. repetitions
is a list of integers denoting the number of subdivisions in each coordinate direction.
If the colorize
flag is set, the boundary_indicators
of the surfaces 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. 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. For instance, the center point (1,-1,1) yields a material id 5.
static 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 | |||
) | [inline, static] |
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 widths 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
.
static void GridGenerator::subdivided_hyper_rectangle | ( | Triangulation< dim > & | tria, | |
const std::vector< std::vector< double > > & | spacing, | |||
const Point< dim > & | p, | |||
const Table< dim, unsigned char > & | material_id, | |||
const bool | colorize = false | |||
) | [inline, static] |
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.
static void GridGenerator::parallelogram | ( | Triangulation< dim > & | tria, | |
const Tensor< 2, dim > & | corners, | |||
const bool | colorize = false | |||
) | [inline, static] |
A parallelogram. The first corner point is the origin. The dim
adjacent points are the one-dimensional subtensors of the tensor provided and additional points will be sums of these two vectors. Colorizing is done according to hyper_rectangle().
static void GridGenerator::enclosed_hyper_cube | ( | Triangulation< dim > & | tria, | |
const double | left = 0. , |
|||
const double | right = 1. , |
|||
const double | thickness = 1. , |
|||
const bool | colorize = false | |||
) | [inline, static] |
Hypercube with a layer of hypercubes around it. The first two parameters 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 id's according to the following scheme: extending over the inner cube in (+/-) x-direction: 1/2. In y-direction 4/8, in z-direction 16/32. The cells at corners and edges (3d) get these values bitwise or'd.
Presently only available in 2d and 3d.
static void GridGenerator::hyper_ball | ( | Triangulation< dim > & | tria, | |
const Point< dim > & | center = Point< dim >() , |
|||
const double | radius = 1. | |||
) | [inline, static] |
Initialize the given triangulation with a hyperball, i.e. a circle or a ball around center
with given radius
.
In order to avoid degenerate cells at the boundaries, the circle is triangulated by five cells, the ball by seven cells. 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.
static void GridGenerator::half_hyper_ball | ( | Triangulation< dim > & | tria, | |
const Point< dim > & | center = Point< dim >() , |
|||
const double | radius = 1. | |||
) | [inline, static] |
This class 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 appropriate boundary class is HalfHyperBallBoundary, or HyperBallBoundary.
static void GridGenerator::cylinder | ( | Triangulation< dim > & | tria, | |
const double | radius = 1. , |
|||
const double | half_length = 1. | |||
) | [inline, static] |
Create a cylinder around the x-axis. 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.
static void GridGenerator::hyper_L | ( | Triangulation< dim > & | tria, | |
const double | left = -1. , |
|||
const double | right = 1. | |||
) | [inline, static] |
Initialize the given triangulation with a hyper-L consisting of exactly 2^dim-1
cells. It produces the hypercube with the interval [left,right] without the hypercube made out of the interval [(a+b)/2,b].
The triangulation needs to be void upon calling this function.
This function is declared to exist for triangulations of all space dimensions, but throws an error if called in 1d.
static void GridGenerator::hyper_cube_slit | ( | Triangulation< dim > & | tria, | |
const double | left = 0. , |
|||
const double | right = 1. , |
|||
const bool | colorize = false | |||
) | [inline, static] |
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.
static void GridGenerator::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 | |||
) | [inline, static] |
Produce a hyper-shell, the region between two spheres around center
, with given inner_radius
and outer_radius
.
If the flag colorize
is true
, then the outer boundary will have the id 1, while the inner boundary has id zero. If the flag is false
, both have id zero.
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 two different numbers are meaningful, 6 for a surface based on a hexahedron and 12 for the rhombic dodecahedron.
This function is declared to exist for triangulations of all space dimensions, but throws an error if called in 1d. It is also currently not implemented in 3d.
static 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 | |||
) | [inline, static] |
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 class 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 is zero (as is the default), then it is computed adaptively such that the resulting elements have the least aspect ratio.
At present, this function only exists in 2d.
static 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 | |||
) | [inline, static] |
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 for this initial triangulation. 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
.
static void GridGenerator::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 | repetition = 1 , |
|||
const bool | colorize = false | |||
) | [inline, static] |
This class produces a square on the xy-plane with a circular hole in the middle, times the interval [0.L] (only in 3d).
It is implemented in 2d and 3d, and takes the following arguments:
inner_radius:
size of the internal hole outer_radius:
size of the biggest enclosed cylinder L:
extension on the z-direction
repetitions:
number of subdivisions along the z-direction
colorize:
wether to assign different boundary indicators to different faces. 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. static void GridGenerator::moebius | ( | Triangulation< 3, 3 > & | tria, | |
const unsigned int | n_cells, | |||
const unsigned int | n_rotations, | |||
const double | R, | |||
const double | r | |||
) | [static] |
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 glueing 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 bend together as loop. |
static void GridGenerator::laplace_transformation | ( | Triangulation< dim > & | tria, | |
const std::map< unsigned int, Point< dim > > & | new_points | |||
) | [inline, static] |
This function transformes the Triangulation
tria
smoothly to a domain that is described by the boundary points in the map new_points
. This map maps the point indices to the boundary points in the transformed domain.
Note, that the Triangulation
is changed in-place, therefore you don't need to keep two triangulations, but the given triangulation is changed (overwritten).
In 1d, this function is not currently implemented.
static void GridGenerator::colorize_hyper_rectangle | ( | Triangulation< dim, spacedim > & | tria | ) | [inline, static, private] |
Perform the action specified by the colorize
flag of the hyper_rectangle() function of this class.
static void GridGenerator::colorize_subdivided_hyper_rectangle | ( | Triangulation< dim > & | tria, | |
const Point< dim > & | p1, | |||
const Point< dim > & | p2, | |||
const double | epsilon | |||
) | [inline, static, private] |
Perform the action specified by the colorize
flag of the subdivided_hyper_rectangle() function of this class. This function is singled out because it is dimension specific.
static void GridGenerator::colorize_hyper_shell | ( | Triangulation< dim > & | tria, | |
const Point< dim > & | center, | |||
const double | inner_radius, | |||
const double | outer_radius | |||
) | [inline, static, private] |
Assign boundary number zero to the inner shell boundary and 1 to the outer.
static void GridGenerator::laplace_solve | ( | const SparseMatrix< double > & | S, | |
const std::map< unsigned int, double > & | m, | |||
Vector< double > & | u | |||
) | [static, private] |
Solve the Laplace equation for laplace_transformation
function for one of the dim
space dimensions. Externalized into a function of its own in order to allow parallel execution.