![]() |
Reference documentation for deal.II version 8.1.0
|
Functions | |
template<int dim, int spacedim> | |
void | hyper_cube (Triangulation< dim, spacedim > &tria, const double left=0., const double right=1.) |
template<int dim> | |
void | subdivided_hyper_cube (Triangulation< dim > &tria, const unsigned int repetitions, const double left=0., const double right=1.) |
template<int dim, int spacedim> | |
void | hyper_rectangle (Triangulation< dim, spacedim > &tria, const Point< spacedim > &p1, const Point< spacedim > &p2, const bool colorize=false) |
template<int dim> | |
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> | |
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> | |
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> | |
void | parallelogram (Triangulation< dim > &tria, const Point< dim > corners[dim], const bool colorize=false) |
template<int dim> | |
void | parallelogram (Triangulation< dim > &tria, const Tensor< 2, dim > &corners, const bool colorize=false) DEAL_II_DEPRECATED |
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> | |
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.) |
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 | 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> | |
void | hyper_L (Triangulation< dim > &tria, const double left=-1., const double right=1.) |
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> | |
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> | |
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) |
void | torus (Triangulation< 2, 3 > &tria, const double R, const double r) |
template<int dim> | |
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) |
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 | merge_triangulations (const Triangulation< dim, spacedim > &triangulation_1, const Triangulation< dim, spacedim > &triangulation_2, Triangulation< dim, spacedim > &result) |
void | extrude_triangulation (const Triangulation< 2, 2 > &input, const unsigned int n_slices, const double height, Triangulation< 3, 3 > &result) |
template<int dim> | |
void | laplace_transformation (Triangulation< dim > &tria, const std::map< unsigned int, Point< dim > > &new_points) DEAL_II_DEPRECATED |
DeclException0 (ExcInvalidRadii) | |
DeclException1 (ExcInvalidRepetitions, int,<< "The number of repetitions "<< arg1<< " must be >=1.") | |
DeclException1 (ExcInvalidRepetitionsDimension, int,<< "The vector of repetitions must have "<< arg1<<" elements.") | |
This namespace provides a collection of functions for generating triangulations for some basic geometries.
Some of these functions receive a flag colorize
. If this is set, parts of the boundary receive different boundary indicators (GlossBoundaryIndicator), allowing them to be distinguished for the purpose of attaching geometry objects and evaluating different boundary conditions.
This namespace also provides a function GridGenerator::laplace_transformation that smoothly transforms a domain into another one. This can be used to transform basic geometries to more complicated ones, like a shell to a grid of an airfoil, for example.
void GridGenerator::hyper_cube | ( | Triangulation< dim, spacedim > & | tria, |
const double | left = 0. , |
||
const double | right = 1. |
||
) |
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. All boundary indicators are set to zero ("not colorized") for 2d and 3d. In 1d the indicators are colorized, see hyper_rectangle().
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.
void GridGenerator::subdivided_hyper_cube | ( | Triangulation< dim > & | tria, |
const unsigned int | repetitions, | ||
const double | left = 0. , |
||
const double | right = 1. |
||
) |
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.
void GridGenerator::hyper_rectangle | ( | Triangulation< dim, spacedim > & | tria, |
const Point< spacedim > & | p1, | ||
const Point< spacedim > & | 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 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.
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 |
||
) |
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 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. 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.
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 | ||
) |
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
.
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.
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 dim
adjacent points are the ones given in the second argument and the fourth point will be the sum of these two vectors. Colorizing is done in the same way as in hyper_rectangle().
void GridGenerator::parallelogram | ( | Triangulation< dim > & | tria, |
const Tensor< 2, dim > & | corners, | ||
const bool | 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().
GridReordering::reorder_grid
). In other words, if reodering of the vertices does occur, the ordering of vertices in the array of corners
will no longer refer to the same triangulation.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().
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, ie. 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().
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. 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.
void GridGenerator::hyper_ball | ( | Triangulation< dim > & | tria, |
const Point< dim > & | center = Point< dim >() , |
||
const double | radius = 1. |
||
) |
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.
void GridGenerator::half_hyper_ball | ( | Triangulation< dim > & | tria, |
const Point< dim > & | center = Point< dim >() , |
||
const double | radius = 1. |
||
) |
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.
void GridGenerator::cylinder | ( | Triangulation< dim > & | tria, |
const double | radius = 1. , |
||
const double | half_length = 1. |
||
) |
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.
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 cutted 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.
An example of use can be found in the documentation of the ConeBoundary class, with which you probably want to associate boundary indicator 0 (the hull of the cone).
void GridGenerator::hyper_L | ( | Triangulation< dim > & | tria, |
const double | left = -1. , |
||
const double | right = 1. |
||
) |
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.
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.
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 |
||
) |
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. If the flag is false
, both have indicator 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 (i.e. 6 panels on the inner sphere extruded in radial direction to form 6 cells) and 12 for the rhombic dodecahedron. These give rise to the following meshes upon one refinement:
Neither of these meshes is particularly good since one ends up with poorly shaped cells at the inner edge upon refinement. For example, this is the middle plane of the mesh for the n_cells=6
:
The mesh generated with n_cells=6
is better but still not good. As a consequence, you may also specify n_cells=96
as a third option. The mesh generated in this way is based on a once refined version of the one with n_cells=12
, where all internal nodes are re-placed along a shell somewhere between the inner and outer boundary of the domain. The following two images compare half of the hyper shell for n_cells=12
and n_cells=96
(note that the doubled radial lines on the cross section are artifacts of the visualization):
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 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.
If colorize is set to true, the inner, outer, left, and right boundary get indicator 0, 1, 2, and 3, respectively. Otherwise all indicators are set to 0.
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.
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 |
||
) |
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
.
void GridGenerator::torus | ( | Triangulation< 2, 3 > & | tria, |
const double | R, | ||
const double | r | ||
) |
Produce the surface meshing of the torus. The axis of the torus is the -axis while the plane of the torus is the
-
plane. The boundary of this object can be described by the TorusBoundary class.
tria | The triangulation to be filled. |
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 inner radius of the torus. |
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 |
||
) |
This class produces a square on the xy-plane with a circular hole in the middle. Square and circle are centered at the origin. In 3d, this geometry is extruded in direction to the interval
.
It is implemented in 2d and 3d, and takes the following arguments:
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. 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. 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 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. |
void GridGenerator::merge_triangulations | ( | 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 cells of both triangulation and store it in the third parameter. Previous content of result
will be deleted.
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).
void GridGenerator::extrude_triangulation | ( | const Triangulation< 2, 2 > & | input, |
const unsigned int | n_slices, | ||
const double | height, | ||
Triangulation< 3, 3 > & | result | ||
) |
Take a 2d Triangulation that is being extruded in z direction by the total height of height
using n_slices
slices (minimum is 2). 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.
void GridGenerator::laplace_transformation | ( | Triangulation< dim > & | tria, |
const std::map< unsigned int, Point< dim > > & | new_points | ||
) |
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.
GridGenerator::DeclException0 | ( | ExcInvalidRadii | ) |
Exception
GridGenerator::DeclException1 | ( | ExcInvalidRepetitions | , |
int | , | ||
<< "The number of repetitions "<< arg1<< " must be >=1." | |||
) |
Exception
GridGenerator::DeclException1 | ( | ExcInvalidRepetitionsDimension | , |
int | , | ||
<< "The vector of repetitions must have "<< arg1<<" elements." | |||
) |
Exception