OGS
MeshToolsLib::MeshGenerator Namespace Reference

Namespaces

namespace  AddFaultToVoxelGrid
namespace  VoxelGridFromMesh

Functions

std::vector< MeshLib::Node * > generateRegularPyramidTopNodes (std::vector< double > const &x_coords, std::vector< double > const &y_coords, std::vector< double > const &z_coords, const MathLib::Point3d &origin)
std::vector< MeshLib::Node * > generateRegularNodes (const std::vector< const std::vector< double > * > &vec_xyz_coords, const MathLib::Point3d &origin=MathLib::ORIGIN)
std::vector< MeshLib::Node * > generateRegularNodes (const std::vector< double > &vec_x_coords, const MathLib::Point3d &origin=MathLib::ORIGIN)
std::vector< MeshLib::Node * > generateRegularNodes (std::vector< double > &vec_x_coords, std::vector< double > &vec_y_coords, const MathLib::Point3d &origin=MathLib::ORIGIN)
std::vector< MeshLib::Node * > generateRegularNodes (std::vector< double > &vec_x_coords, std::vector< double > &vec_y_coords, std::vector< double > &vec_z_coords, const MathLib::Point3d &origin=MathLib::ORIGIN)
std::vector< MeshLib::Node * > generateRegularNodes (const std::array< unsigned, 3 > &n_cells, const std::array< double, 3 > &cell_size, const MathLib::Point3d &origin)
MeshLib::MeshgenerateLineMesh (const BaseLib::ISubdivision &div, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateLineMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateLineMesh (const unsigned n_cells, const double cell_size, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularQuadMesh (const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularQuadMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularQuadMesh (const double x_length, const double y_length, const std::size_t x_subdivision, const std::size_t y_subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularQuadMesh (const unsigned n_x_cells, const unsigned n_y_cells, const double cell_size, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularQuadMesh (const unsigned n_x_cells, const unsigned n_y_cells, const double cell_size_x, const double cell_size_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularHexMesh (const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, const BaseLib::ISubdivision &div_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularHexMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularHexMesh (const double x_length, const double y_length, const double z_length, const std::size_t x_subdivision, const std::size_t y_subdivision, const std::size_t z_subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularHexMesh (const unsigned n_x_cells, const unsigned n_y_cells, const unsigned n_z_cells, const double cell_size, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularHexMesh (const unsigned n_x_cells, const unsigned n_y_cells, const unsigned n_z_cells, const double cell_size_x, const double cell_size_y, const double cell_size_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularPyramidMesh (const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, const BaseLib::ISubdivision &div_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTriMesh (const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTriMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTriMesh (const double x_length, const double y_length, const std::size_t x_subdivision, const std::size_t y_subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTriMesh (const unsigned n_x_cells, const unsigned n_y_cells, const double cell_size, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTriMesh (const unsigned n_x_cells, const unsigned n_y_cells, const double cell_size_x, const double cell_size_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularPrismMesh (const double x_length, const double y_length, const double z_length, const std::size_t x_subdivision, const std::size_t y_subdivision, const std::size_t z_subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularPrismMesh (const unsigned n_x_cells, const unsigned n_y_cells, const unsigned n_z_cells, const double cell_size, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularPrismMesh (const unsigned n_x_cells, const unsigned n_y_cells, const unsigned n_z_cells, const double cell_size_x, const double cell_size_y, const double cell_size_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTetMesh (const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, const BaseLib::ISubdivision &div_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTetMesh (const double x_length, const double y_length, const double z_length, const std::size_t x_subdivision, const std::size_t y_subdivision, const std::size_t z_subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshgenerateRegularTetMesh (const unsigned n_x_cells, const unsigned n_y_cells, const unsigned n_z_cells, const double cell_size_x, const double cell_size_y, const double cell_size_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::MeshcreateSurfaceMesh (std::string const &mesh_name, MathLib::Point3d const &ll, MathLib::Point3d const &ur, std::array< std::size_t, 2 > const &n_steps, const std::function< double(double, double)> &f)

Function Documentation

◆ createSurfaceMesh()

MeshLib::Mesh * MeshToolsLib::MeshGenerator::createSurfaceMesh ( std::string const & mesh_name,
MathLib::Point3d const & ll,
MathLib::Point3d const & ur,
std::array< std::size_t, 2 > const & n_steps,
const std::function< double(double, double)> & f )

Constructs a surface mesh approximating a surface in the 3d space given by a function. The surface within the xy-domain \([ll[0], ur[0]] \times [ll[1], ur[1]\) is described using the function \(f(x,y)\).

Definition at line 741 of file MeshGenerator.cpp.

745{
746 std::array<double, 2> const step_size{{(ur[0] - ll[0]) / (n_steps[0] - 1),
747 (ur[1] - ll[1]) / (n_steps[1] - 1)}};
748
749 std::vector<MeshLib::Node*> nodes;
750 for (std::size_t j(0); j < n_steps[1]; ++j)
751 {
752 for (std::size_t i(0); i < n_steps[0]; ++i)
753 {
754 std::size_t const id = i + j * n_steps[1];
755 double const x = ll[0] + i * step_size[0];
756 double const y = ll[1] + j * step_size[1];
757
758 nodes.push_back(new MeshLib::Node({x, y, f(x, y)}, id));
759 }
760 }
761
762 std::vector<MeshLib::Element*> sfc_eles;
763 for (std::size_t j(0); j < n_steps[1] - 1; ++j)
764 {
765 for (std::size_t i(0); i < n_steps[0] - 1; ++i)
766 {
767 std::size_t id_ll(i + j * n_steps[0]);
768 std::size_t id_lr(i + 1 + j * n_steps[0]);
769 std::size_t id_ul(i + (j + 1) * n_steps[0]);
770 std::size_t id_ur(i + 1 + (j + 1) * n_steps[0]);
771 sfc_eles.push_back(
772 new MeshLib::Tri({nodes[id_ll], nodes[id_lr], nodes[id_ur]}));
773 sfc_eles.push_back(
774 new MeshLib::Tri({nodes[id_ll], nodes[id_ur], nodes[id_ul]}));
775 }
776 }
777
778 return new MeshLib::Mesh(mesh_name, nodes, sfc_eles,
779 true /* compute_element_neighbors */);
780}
TemplateElement< MeshLib::TriRule3 > Tri
Definition Tri.h:15

◆ generateLineMesh() [1/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateLineMesh ( const BaseLib::ISubdivision & div,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate an 1D Line-Element mesh. The mesh is generated in x-direction.

Parameters
divSubdivision operator
originOptional mesh's origin (the left-most point) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 172 of file MeshGenerator.cpp.

175{
176 const std::vector<double> vec_x(div());
177 std::vector<MeshLib::Node*> nodes(generateRegularNodes(vec_x, origin));
178
179 // elements
180 const std::size_t n_cells = nodes.size() - 1;
181 std::vector<MeshLib::Element*> elements;
182 elements.reserve(n_cells);
183
184 for (std::size_t i = 0; i < n_cells; i++)
185 {
186 elements.push_back(new MeshLib::Line({nodes[i], nodes[i + 1]}));
187 }
188
189 return new MeshLib::Mesh(mesh_name, nodes, elements,
190 true /* compute_element_neighbors */);
191}
TemplateElement< MeshLib::LineRule2 > Line
Definition Line.h:14
std::vector< MeshLib::Node * > generateRegularNodes(const std::vector< const std::vector< double > * > &vec_xyz_coords, const MathLib::Point3d &origin=MathLib::ORIGIN)

References generateRegularNodes().

Referenced by CreateStructuredGridDialog::accept(), generateLineMesh(), generateLineMesh(), and main().

◆ generateLineMesh() [2/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateLineMesh ( const double length,
const std::size_t subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate an 1D Line-Element mesh. The mesh is generated in x-direction.

Parameters
lengthMesh's length in x-direction.
subdivisionNumber of subdivisions.
originOptional mesh's origin (the left-most point) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 153 of file MeshGenerator.cpp.

157{
158 return generateLineMesh(subdivision, length / subdivision, origin,
159 mesh_name);
160}
MeshLib::Mesh * generateLineMesh(const BaseLib::ISubdivision &div, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")

References generateLineMesh().

◆ generateLineMesh() [3/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateLineMesh ( const unsigned n_cells,
const double cell_size,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate an 1D Line-Element mesh. The mesh is generated in x-direction.

Parameters
n_cellsNumber of cells.
cell_sizeLength of Line elements
originOptional mesh's origin (the left-most point) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 162 of file MeshGenerator.cpp.

166{
167 return generateLineMesh(
168 BaseLib::UniformSubdivision(n_cells * cell_size, n_cells), origin,
169 mesh_name);
170}

References generateLineMesh().

◆ generateRegularHexMesh() [1/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularHexMesh ( const BaseLib::ISubdivision & div_x,
const BaseLib::ISubdivision & div_y,
const BaseLib::ISubdivision & div_z,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Hex-Element mesh.

Parameters
div_xSubdivision operator in x direction
div_ySubdivision operator in y direction
div_zSubdivision operator in z direction
originOptional mesh's origin (the bottom lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 332 of file MeshGenerator.cpp.

338{
339 std::vector<double> vec_x(div_x());
340 std::vector<double> vec_y(div_y());
341 std::vector<double> vec_z(div_z());
342 std::vector<MeshLib::Node*> nodes(
343 generateRegularNodes(vec_x, vec_y, vec_z, origin));
344
345 const unsigned n_x_nodes(vec_x.size());
346 const unsigned n_y_nodes(vec_y.size());
347 const unsigned n_x_cells(vec_x.size() - 1);
348 const unsigned n_y_cells(vec_y.size() - 1);
349 const unsigned n_z_cells(vec_z.size() - 1);
350
351 // elements
352 std::vector<MeshLib::Element*> elements;
353 elements.reserve(n_x_cells * n_y_cells * n_z_cells);
354
355 for (std::size_t i = 0; i < n_z_cells; i++)
356 {
357 const std::size_t offset_z1 = i * n_x_nodes * n_y_nodes; // bottom
358 const std::size_t offset_z2 = (i + 1) * n_x_nodes * n_y_nodes; // top
359 for (std::size_t j = 0; j < n_y_cells; j++)
360 {
361 const std::size_t offset_y1 = j * n_x_nodes;
362 const std::size_t offset_y2 = (j + 1) * n_x_nodes;
363 for (std::size_t k = 0; k < n_x_cells; k++)
364 {
365 elements.push_back(
366 new MeshLib::Hex({// bottom
367 nodes[offset_z1 + offset_y1 + k],
368 nodes[offset_z1 + offset_y1 + k + 1],
369 nodes[offset_z1 + offset_y2 + k + 1],
370 nodes[offset_z1 + offset_y2 + k],
371 // top
372 nodes[offset_z2 + offset_y1 + k],
373 nodes[offset_z2 + offset_y1 + k + 1],
374 nodes[offset_z2 + offset_y2 + k + 1],
375 nodes[offset_z2 + offset_y2 + k]}));
376 }
377 }
378 }
379
380 return new MeshLib::Mesh(mesh_name, nodes, elements,
381 true /* compute_element_neighbors */);
382}
TemplateElement< MeshLib::HexRule8 > Hex
Definition Hex.h:14

References generateRegularNodes().

Referenced by CreateStructuredGridDialog::accept(), Vtu2GridDialog::accept(), MeshToolsLib::RasterToMesh::convert(), generateInitialMesh(), generateRegularHexMesh(), generateRegularHexMesh(), generateRegularHexMesh(), generateRegularHexMesh(), and main().

◆ generateRegularHexMesh() [2/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularHexMesh ( const double length,
const std::size_t subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Hex-Element mesh.

Parameters
lengthMesh dimensions in x- and y- and z-directions.
subdivisionNumber of subdivisions.
originOptional mesh's origin (the bottom lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 276 of file MeshGenerator.cpp.

281{
283 subdivision, subdivision, subdivision, length / subdivision, origin,
284 mesh_name);
285}
MeshLib::Mesh * generateRegularHexMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, const BaseLib::ISubdivision &div_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")

References generateRegularHexMesh().

◆ generateRegularHexMesh() [3/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularHexMesh ( const double x_length,
const double y_length,
const double z_length,
const std::size_t x_subdivision,
const std::size_t y_subdivision,
const std::size_t z_subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Hex-Element mesh.

Parameters
x_lengthMesh dimension in x-direction.
y_lengthMesh dimension in y-direction.
z_lengthMesh dimension in z-direction.
x_subdivisionNumber of subdivisions in x-direction.
y_subdivisionNumber of subdivisions in y-direction.
z_subdivisionNumber of subdivisions in z-direction.
originOptional mesh's origin (the bottom lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 287 of file MeshGenerator.cpp.

296{
298 x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
299 y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
300}

References generateRegularHexMesh().

◆ generateRegularHexMesh() [4/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularHexMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const unsigned n_z_cells,
const double cell_size,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Hex-Element mesh.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
n_z_cellsNumber of cells in z-direction.
cell_sizeEdge length of Hex elements
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 302 of file MeshGenerator.cpp.

309{
311 n_x_cells, n_y_cells, n_z_cells, cell_size, cell_size, cell_size,
312 origin, mesh_name);
313}

References generateRegularHexMesh().

◆ generateRegularHexMesh() [5/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularHexMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const unsigned n_z_cells,
const double cell_size_x,
const double cell_size_y,
const double cell_size_z,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Hex-Element mesh.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
n_z_cellsNumber of cells in z-direction.
cell_size_xEdge length of Hex elements in x-direction.
cell_size_yEdge length of Hex elements in y s-direction.
cell_size_zEdge length of Hex elements in z-direction
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 315 of file MeshGenerator.cpp.

324{
326 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
327 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells),
328 BaseLib::UniformSubdivision(n_z_cells * cell_size_z, n_z_cells), origin,
329 mesh_name);
330}

References generateRegularHexMesh().

◆ generateRegularNodes() [1/5]

std::vector< MeshLib::Node * > MeshToolsLib::MeshGenerator::generateRegularNodes ( const std::array< unsigned, 3 > & n_cells,
const std::array< double, 3 > & cell_size,
const MathLib::Point3d & origin )

Generate regularly placed mesh nodes in 3D spaces

Parameters
n_cellsan array of the number of cells in x,y,z directions
cell_sizean array of cell sizes in x,y,z directions
origincoordinates of the left-most point
Returns
a vector of created mesh nodes

Definition at line 96 of file MeshGenerator.cpp.

100{
101 std::vector<MeshLib::Node*> nodes;
102 nodes.reserve((n_cells[0] + 1) * (n_cells[1] + 1) * (n_cells[2] + 1));
103
104 for (std::size_t i = 0; i < n_cells[2] + 1; i++)
105 {
106 const double z(origin[2] + cell_size[2] * i);
107 for (std::size_t j = 0; j < n_cells[1] + 1; j++)
108 {
109 const double y(origin[1] + cell_size[1] * j);
110 for (std::size_t k = 0; k < n_cells[0] + 1; k++)
111 {
112 nodes.push_back(
113 new MeshLib::Node(origin[0] + cell_size[0] * k, y, z));
114 }
115 }
116 }
117 return nodes;
118}

◆ generateRegularNodes() [2/5]

std::vector< MeshLib::Node * > MeshToolsLib::MeshGenerator::generateRegularNodes ( const std::vector< const std::vector< double > * > & vec_xyz_coords,
const MathLib::Point3d & origin = MathLib::ORIGIN )

Generate regularly placed mesh nodes in 3D spaces

Parameters
vec_xyz_coordsa vector of coordinates in x,y,z directions
origincoordinates of the left-most point
Returns
a vector of created mesh nodes

Definition at line 21 of file MeshGenerator.cpp.

24{
25 auto const shift_coordinates =
26 [](auto const& in, auto& out, auto const& shift)
27 {
28 std::transform(in.begin(), in.end(), std::back_inserter(out),
29 [&shift](auto const& v) { return v + shift; });
30 };
31 std::array<std::vector<double>, 3> coords;
32 for (std::size_t i = 0; i < 3; ++i)
33 {
34 coords[i].reserve(vec_xyz_coords[i]->size());
35 shift_coordinates(*vec_xyz_coords[i], coords[i], origin[i]);
36 }
37
38 std::vector<MeshLib::Node*> nodes;
39 nodes.reserve(coords[0].size() * coords[1].size() * coords[2].size());
40
41 for (auto const z : coords[2])
42 {
43 for (auto const y : coords[1])
44 {
45 std::transform(coords[0].begin(), coords[0].end(),
46 std::back_inserter(nodes),
47 [&y, &z](double const& x)
48 { return new MeshLib::Node(x, y, z); });
49 }
50 }
51 return nodes;
52}

References MeshLib::Node.

Referenced by generateLineMesh(), generateRegularHexMesh(), generateRegularNodes(), generateRegularNodes(), generateRegularNodes(), generateRegularPyramidMesh(), generateRegularQuadMesh(), generateRegularTetMesh(), and generateRegularTriMesh().

◆ generateRegularNodes() [3/5]

std::vector< MeshLib::Node * > MeshToolsLib::MeshGenerator::generateRegularNodes ( const std::vector< double > & vec_x_coords,
const MathLib::Point3d & origin = MathLib::ORIGIN )

Generate regularly placed mesh nodes in 1D space

Parameters
vec_x_coordsa vector of x coordinates
origincoordinates of the left-most point
Returns
a vector of created mesh nodes

Definition at line 54 of file MeshGenerator.cpp.

56{
57 std::vector<const std::vector<double>*> vec_xyz_coords;
58 vec_xyz_coords.push_back(&vec_x_coords);
59 std::vector<double> dummy(1, 0.0);
60 for (unsigned i = vec_xyz_coords.size() - 1; i < 3u; i++)
61 {
62 vec_xyz_coords.push_back(&dummy);
63 }
64 return generateRegularNodes(vec_xyz_coords, origin);
65}

References generateRegularNodes().

◆ generateRegularNodes() [4/5]

std::vector< MeshLib::Node * > MeshToolsLib::MeshGenerator::generateRegularNodes ( std::vector< double > & vec_x_coords,
std::vector< double > & vec_y_coords,
const MathLib::Point3d & origin = MathLib::ORIGIN )

Generate regularly placed mesh nodes in 1D space

Parameters
vec_x_coordsa vector of x coordinates
vec_y_coordsa vector of y coordinates
origincoordinates of the left-most point
Returns
a vector of created mesh nodes

Definition at line 67 of file MeshGenerator.cpp.

71{
72 std::vector<const std::vector<double>*> vec_xyz_coords;
73 vec_xyz_coords.push_back(&vec_x_coords);
74 vec_xyz_coords.push_back(&vec_y_coords);
75 std::vector<double> dummy(1, 0.0);
76 for (unsigned i = vec_xyz_coords.size() - 1; i < 3u; i++)
77 {
78 vec_xyz_coords.push_back(&dummy);
79 }
80 return generateRegularNodes(vec_xyz_coords, origin);
81}

References generateRegularNodes().

◆ generateRegularNodes() [5/5]

std::vector< MeshLib::Node * > MeshToolsLib::MeshGenerator::generateRegularNodes ( std::vector< double > & vec_x_coords,
std::vector< double > & vec_y_coords,
std::vector< double > & vec_z_coords,
const MathLib::Point3d & origin = MathLib::ORIGIN )

Generate regularly placed mesh nodes in 1D space

Parameters
vec_x_coordsa vector of x coordinates
vec_y_coordsa vector of y coordinates
vec_z_coordsa vector of z coordinates
origincoordinates of the left-most point
Returns
a vector of created mesh nodes

Definition at line 83 of file MeshGenerator.cpp.

88{
89 std::vector<const std::vector<double>*> vec_xyz_coords;
90 vec_xyz_coords.push_back(&vec_x_coords);
91 vec_xyz_coords.push_back(&vec_y_coords);
92 vec_xyz_coords.push_back(&vec_z_coords);
93 return generateRegularNodes(vec_xyz_coords, origin);
94}

References generateRegularNodes().

◆ generateRegularPrismMesh() [1/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularPrismMesh ( const double x_length,
const double y_length,
const double z_length,
const std::size_t x_subdivision,
const std::size_t y_subdivision,
const std::size_t z_subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Prism-Element mesh.

Parameters
x_lengthMesh's dimension in x-direction.
y_lengthMesh's dimension in y-direction.
z_lengthMesh's dimension in z-direction.
x_subdivisionNumber of subdivisions in x-direction.
y_subdivisionNumber of subdivisions in y-direction.
z_subdivisionNumber of subdivisions in z-direction.
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 568 of file MeshGenerator.cpp.

577{
579 x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
580 y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
581}
MeshLib::Mesh * generateRegularPrismMesh(const double x_length, const double y_length, const double z_length, const std::size_t x_subdivision, const std::size_t y_subdivision, const std::size_t z_subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")

References generateRegularPrismMesh().

Referenced by CreateStructuredGridDialog::accept(), MeshToolsLib::RasterToMesh::convert(), generateRegularPrismMesh(), generateRegularPrismMesh(), and main().

◆ generateRegularPrismMesh() [2/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularPrismMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const unsigned n_z_cells,
const double cell_size,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Prism-Element mesh.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
n_z_cellsNumber of cells in z-direction.
cell_sizeEdge length of two equal sides of isosceles triangles + height of prism
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 583 of file MeshGenerator.cpp.

590{
591 return generateRegularPrismMesh(n_x_cells, n_y_cells, n_z_cells, cell_size,
592 cell_size, cell_size, origin, mesh_name);
593}

References generateRegularPrismMesh().

◆ generateRegularPrismMesh() [3/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularPrismMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const unsigned n_z_cells,
const double cell_size_x,
const double cell_size_y,
const double cell_size_z,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Prism-Element mesh.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
n_z_cellsNumber of cells in z-direction.
cell_size_xEdge length of triangles in x-direction.
cell_size_yEdge length of triangles in y-direction.
cell_size_zEdge length of triangles in z-direction.
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 595 of file MeshGenerator.cpp.

604{
605 std::unique_ptr<MeshLib::Mesh> mesh(generateRegularTriMesh(
606 n_x_cells, n_y_cells, cell_size_x, cell_size_y, origin, mesh_name));
607 std::size_t const n_tris(mesh->getNumberOfElements());
608 bool const copy_material_ids = false;
609 for (std::size_t i = 0; i < n_z_cells; ++i)
610 {
611 mesh.reset(MeshToolsLib::addLayerToMesh(*mesh, cell_size_z, mesh_name,
612 true, copy_material_ids,
613 std::nullopt));
614 }
615 std::vector<std::size_t> elem_ids(n_tris);
616 std::iota(elem_ids.begin(), elem_ids.end(), 0);
617 return MeshToolsLib::removeElements(*mesh, elem_ids, mesh_name);
618}
MeshLib::Mesh * generateRegularTriMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
MeshLib::Mesh * addLayerToMesh(MeshLib::Mesh const &mesh, double const thickness, std::string const &name, bool const on_top, bool const copy_material_ids, std::optional< int > const layer_id)
MeshLib::Mesh * removeElements(const MeshLib::Mesh &mesh, const std::vector< std::size_t > &removed_element_ids, const std::string &new_mesh_name)

References MeshToolsLib::addLayerToMesh(), generateRegularTriMesh(), and MeshToolsLib::removeElements().

◆ generateRegularPyramidMesh()

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularPyramidMesh ( const BaseLib::ISubdivision & div_x,
const BaseLib::ISubdivision & div_y,
const BaseLib::ISubdivision & div_z,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D mesh that consists of pyramid elements.

This algorithm generates regular grid points as well as middle points of the virtual hexahedra cells and split each hex grid cell into six pyramids.

Parameters
div_xSubdivision operator in x direction
div_ySubdivision operator in y direction
div_zSubdivision operator in z direction
originOptional mesh's origin (the bottom lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 384 of file MeshGenerator.cpp.

390{
391 std::vector<double> vec_x(div_x());
392 std::vector<double> vec_y(div_y());
393 std::vector<double> vec_z(div_z());
394 std::vector<MeshLib::Node*> nodes(
395 generateRegularNodes(vec_x, vec_y, vec_z, origin));
396 std::vector<MeshLib::Node*> const top_nodes(
397 generateRegularPyramidTopNodes(vec_x, vec_y, vec_z, origin));
398
399 nodes.insert(nodes.end(), top_nodes.begin(), top_nodes.end());
400
401 const unsigned n_x_nodes(vec_x.size());
402 const unsigned n_y_nodes(vec_y.size());
403 const unsigned n_z_nodes(vec_z.size());
404 const unsigned n_x_cells(vec_x.size() - 1);
405 const unsigned n_y_cells(vec_y.size() - 1);
406 const unsigned n_z_cells(vec_z.size() - 1);
407
408 // elements
409 std::vector<MeshLib::Element*> elements;
410 auto const top_node_offset(n_x_nodes * n_y_nodes * n_z_nodes);
411 elements.reserve(n_x_cells * n_y_cells * n_z_cells);
412
413 for (std::size_t i = 0; i < n_z_cells; i++)
414 {
415 const std::size_t offset_z1 = i * n_x_nodes * n_y_nodes; // bottom
416 const std::size_t offset_z2 = (i + 1) * n_x_nodes * n_y_nodes; // top
417 for (std::size_t j = 0; j < n_y_cells; j++)
418 {
419 const std::size_t offset_y1 = j * n_x_nodes;
420 const std::size_t offset_y2 = (j + 1) * n_x_nodes;
421 for (std::size_t k = 0; k < n_x_cells; k++)
422 {
423 // generate 6 pyramids within the virtual hexahedron cell
424 int const pyramid_top_index(i * n_x_cells * n_y_cells +
425 j * n_x_cells + k +
426 top_node_offset);
427 elements.push_back(
428 new MeshLib::Pyramid{{// bottom 'hexahedron' face
429 nodes[offset_z1 + offset_y1 + k],
430 nodes[offset_z1 + offset_y1 + k + 1],
431 nodes[offset_z1 + offset_y2 + k + 1],
432 nodes[offset_z1 + offset_y2 + k],
433 // top
434 nodes[pyramid_top_index]}});
435 elements.push_back(new MeshLib::Pyramid{
436 {// top 'hexahedron' face
437 nodes[offset_z2 + offset_y1 + k + 1],
438 nodes[offset_z2 + offset_y1 + k],
439 nodes[offset_z2 + offset_y2 + k],
440 nodes[offset_z2 + offset_y2 + k + 1],
441 // top of pyramid directed towards the bottom
442 nodes[pyramid_top_index]}});
443 elements.push_back(new MeshLib::Pyramid{
444 {// right 'hexahedron' face
445 nodes[offset_z1 + offset_y1 + k + 1],
446 nodes[offset_z2 + offset_y1 + k + 1],
447 nodes[offset_z2 + offset_y2 + k + 1],
448 nodes[offset_z1 + offset_y2 + k + 1],
449 // top of pyramid directed towards the bottom
450 nodes[pyramid_top_index]}});
451 elements.push_back(new MeshLib::Pyramid{
452 {// left 'hexahedron' face
453 nodes[offset_z2 + offset_y1 + k],
454 nodes[offset_z1 + offset_y1 + k],
455 nodes[offset_z1 + offset_y2 + k],
456 nodes[offset_z2 + offset_y2 + k],
457 // top of pyramid directed towards the bottom
458 nodes[pyramid_top_index]}});
459 elements.push_back(new MeshLib::Pyramid{
460 {// front 'hexahedron' face
461 nodes[offset_z2 + offset_y1 + k],
462 nodes[offset_z2 + offset_y1 + k + 1],
463 nodes[offset_z1 + offset_y1 + k + 1],
464 nodes[offset_z1 + offset_y1 + k],
465 // top of pyramid directed towards the bottom
466 nodes[pyramid_top_index]}});
467 elements.push_back(new MeshLib::Pyramid{
468 {// back 'hexahedron' face
469 nodes[offset_z1 + offset_y2 + k],
470 nodes[offset_z1 + offset_y2 + k + 1],
471 nodes[offset_z2 + offset_y2 + k + 1],
472 nodes[offset_z2 + offset_y2 + k],
473 // top of pyramid directed towards the bottom
474 nodes[pyramid_top_index]}});
475 }
476 }
477 }
478 return new MeshLib::Mesh(mesh_name, nodes, elements,
479 true /* compute_element_neighbors */);
480}
TemplateElement< MeshLib::PyramidRule5 > Pyramid
Definition Pyramid.h:14
std::vector< MeshLib::Node * > generateRegularPyramidTopNodes(std::vector< double > const &x_coords, std::vector< double > const &y_coords, std::vector< double > const &z_coords, const MathLib::Point3d &origin)

References generateRegularNodes(), and generateRegularPyramidTopNodes().

Referenced by main().

◆ generateRegularPyramidTopNodes()

std::vector< MeshLib::Node * > MeshToolsLib::MeshGenerator::generateRegularPyramidTopNodes ( std::vector< double > const & x_coords,
std::vector< double > const & y_coords,
std::vector< double > const & z_coords,
const MathLib::Point3d & origin )

Definition at line 122 of file MeshGenerator.cpp.

127{
128 std::vector<MeshLib::Node*> nodes;
129 nodes.reserve((x_coords.size() - 1) * (y_coords.size() - 1) *
130 (z_coords.size() - 1));
131
132 auto const n_x_coords = x_coords.size() - 1;
133 auto const n_y_coords = y_coords.size() - 1;
134 auto const n_z_coords = z_coords.size() - 1;
135
136 for (std::size_t i = 0; i < n_z_coords; i++)
137 {
138 const double z((z_coords[i] + z_coords[i + 1]) / 2 + origin[2]);
139 for (std::size_t j = 0; j < n_y_coords; j++)
140 {
141 const double y((y_coords[j] + y_coords[j + 1]) / 2 + origin[1]);
142 for (std::size_t k = 0; k < n_x_coords; k++)
143 {
144 const double x((x_coords[k] + x_coords[k + 1]) / 2 + origin[0]);
145 nodes.push_back(new MeshLib::Node(x, y, z));
146 }
147 }
148 }
149 return nodes;
150}

Referenced by generateRegularPyramidMesh().

◆ generateRegularQuadMesh() [1/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularQuadMesh ( const BaseLib::ISubdivision & div_x,
const BaseLib::ISubdivision & div_y,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Quad-Element mesh. The mesh is generated in the x-y-plane.

Parameters
div_xSubdivision operator in x direction
div_ySubdivision operator in y direction
originOptional mesh's origin (the left-most point) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 242 of file MeshGenerator.cpp.

247{
248 std::vector<double> vec_x(div_x());
249 std::vector<double> vec_y(div_y());
250 std::vector<MeshLib::Node*> nodes(
251 generateRegularNodes(vec_x, vec_y, origin));
252 const unsigned n_x_nodes(vec_x.size());
253
254 // elements
255 std::vector<MeshLib::Element*> elements;
256 const unsigned n_x_cells(vec_x.size() - 1);
257 const unsigned n_y_cells(vec_y.size() - 1);
258 elements.reserve(n_x_cells * n_y_cells);
259
260 for (std::size_t j = 0; j < n_y_cells; j++)
261 {
262 const std::size_t offset_y1 = j * n_x_nodes;
263 const std::size_t offset_y2 = (j + 1) * n_x_nodes;
264 for (std::size_t k = 0; k < n_x_cells; k++)
265 {
266 elements.push_back(new MeshLib::Quad(
267 {nodes[offset_y1 + k], nodes[offset_y1 + k + 1],
268 nodes[offset_y2 + k + 1], nodes[offset_y2 + k]}));
269 }
270 }
271
272 return new MeshLib::Mesh(mesh_name, nodes, elements,
273 true /* compute_element_neighbors */);
274}
TemplateElement< MeshLib::QuadRule4 > Quad
Definition Quad.h:17

References generateRegularNodes().

Referenced by CreateStructuredGridDialog::accept(), MeshToolsLib::RasterToMesh::convert(), generateRegularQuadMesh(), generateRegularQuadMesh(), generateRegularQuadMesh(), generateRegularQuadMesh(), and main().

◆ generateRegularQuadMesh() [2/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularQuadMesh ( const double length,
const std::size_t subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Quad-Element mesh. The mesh is generated in the x-y-plane.

Parameters
lengthMesh's dimensions in x- and y-directions.
subdivisionNumber of subdivisions.
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 193 of file MeshGenerator.cpp.

198{
199 return generateRegularQuadMesh(subdivision, subdivision,
200 length / subdivision, length / subdivision,
201 origin, mesh_name);
202}
MeshLib::Mesh * generateRegularQuadMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")

References generateRegularQuadMesh().

◆ generateRegularQuadMesh() [3/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularQuadMesh ( const double x_length,
const double y_length,
const std::size_t x_subdivision,
const std::size_t y_subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Quad-Element mesh. The mesh is generated in the x-y-plane.

Parameters
x_lengthMesh's dimension in x-direction.
y_lengthMesh's dimension in y-direction.
x_subdivisionNumber of subdivisions in x-direction.
y_subdivisionNumber of subdivisions in y-direction.
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 204 of file MeshGenerator.cpp.

211{
212 return generateRegularQuadMesh(x_subdivision, y_subdivision,
213 x_length / x_subdivision,
214 y_length / y_subdivision, origin, mesh_name);
215}

References generateRegularQuadMesh().

◆ generateRegularQuadMesh() [4/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularQuadMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const double cell_size,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Quad-Element mesh. The mesh is generated in the x-y-plane.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
cell_sizeEdge length of Quad elements
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 217 of file MeshGenerator.cpp.

223{
224 return generateRegularQuadMesh(n_x_cells, n_y_cells, cell_size, cell_size,
225 origin, mesh_name);
226}

References generateRegularQuadMesh().

◆ generateRegularQuadMesh() [5/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularQuadMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const double cell_size_x,
const double cell_size_y,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Quad-Element mesh. The mesh is generated in the x-y-plane.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
cell_size_xEdge length of Quad elements in x-direction
cell_size_yEdge length of Quad elements in y-direction
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 228 of file MeshGenerator.cpp.

235{
237 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
238 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells), origin,
239 mesh_name);
240}

References generateRegularQuadMesh().

◆ generateRegularTetMesh() [1/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTetMesh ( const BaseLib::ISubdivision & div_x,
const BaseLib::ISubdivision & div_y,
const BaseLib::ISubdivision & div_z,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Tet-Element mesh.

This algorithm generates regular grid points and split each grid cell into six tetrahedrals.

Parameters
div_xSubdivision operator in x direction
div_ySubdivision operator in y direction
div_zSubdivision operator in z direction
originOptional mesh's origin (the bottom lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 652 of file MeshGenerator.cpp.

658{
659 std::vector<double> vec_x(div_x());
660 std::vector<double> vec_y(div_y());
661 std::vector<double> vec_z(div_z());
662 std::vector<MeshLib::Node*> nodes(
663 generateRegularNodes(vec_x, vec_y, vec_z, origin));
664
665 const unsigned n_x_nodes(vec_x.size());
666 const unsigned n_y_nodes(vec_y.size());
667 const unsigned n_x_cells(vec_x.size() - 1);
668 const unsigned n_y_cells(vec_y.size() - 1);
669 const unsigned n_z_cells(vec_z.size() - 1);
670
671 // elements
672 std::vector<MeshLib::Element*> elements;
673 elements.reserve(n_x_cells * n_y_cells * n_z_cells * 6);
674
675 for (std::size_t i = 0; i < n_z_cells; i++)
676 {
677 const std::size_t offset_z1 = i * n_x_nodes * n_y_nodes; // bottom
678 const std::size_t offset_z2 = (i + 1) * n_x_nodes * n_y_nodes; // top
679 for (std::size_t j = 0; j < n_y_cells; j++)
680 {
681 const std::size_t offset_y1 = j * n_x_nodes;
682 const std::size_t offset_y2 = (j + 1) * n_x_nodes;
683 for (std::size_t k = 0; k < n_x_cells; k++)
684 {
685 // tet 1
686 elements.push_back(
687 new MeshLib::Tet({// bottom
688 nodes[offset_z1 + offset_y1 + k],
689 nodes[offset_z1 + offset_y2 + k + 1],
690 nodes[offset_z1 + offset_y2 + k],
691 // top
692 nodes[offset_z2 + offset_y1 + k]}));
693 // tet 2
694 elements.push_back(
695 new MeshLib::Tet({// bottom
696 nodes[offset_z1 + offset_y2 + k + 1],
697 nodes[offset_z1 + offset_y2 + k],
698 // top
699 nodes[offset_z2 + offset_y1 + k],
700 nodes[offset_z2 + offset_y2 + k + 1]}));
701 // tet 3
702 elements.push_back(
703 new MeshLib::Tet({// bottom
704 nodes[offset_z1 + offset_y2 + k],
705 // top
706 nodes[offset_z2 + offset_y1 + k],
707 nodes[offset_z2 + offset_y2 + k + 1],
708 nodes[offset_z2 + offset_y2 + k]}));
709 // tet 4
710 elements.push_back(
711 new MeshLib::Tet({// bottom
712 nodes[offset_z1 + offset_y1 + k],
713 nodes[offset_z1 + offset_y1 + k + 1],
714 nodes[offset_z1 + offset_y2 + k + 1],
715 // top
716 nodes[offset_z2 + offset_y1 + k + 1]}));
717 // tet 5
718 elements.push_back(
719 new MeshLib::Tet({// bottom
720 nodes[offset_z1 + offset_y1 + k],
721 nodes[offset_z1 + offset_y2 + k + 1],
722 // top
723 nodes[offset_z2 + offset_y1 + k],
724 nodes[offset_z2 + offset_y1 + k + 1]}));
725 // tet 6
726 elements.push_back(
727 new MeshLib::Tet({// bottom
728 nodes[offset_z1 + offset_y2 + k + 1],
729 // top
730 nodes[offset_z2 + offset_y1 + k],
731 nodes[offset_z2 + offset_y1 + k + 1],
732 nodes[offset_z2 + offset_y2 + k + 1]}));
733 }
734 }
735 }
736
737 return new MeshLib::Mesh(mesh_name, nodes, elements,
738 true /* compute_element_neighbors */);
739}
TemplateElement< MeshLib::TetRule4 > Tet
Definition Tet.h:14

References generateRegularNodes().

Referenced by generateRegularTetMesh(), generateRegularTetMesh(), and main().

◆ generateRegularTetMesh() [2/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTetMesh ( const double x_length,
const double y_length,
const double z_length,
const std::size_t x_subdivision,
const std::size_t y_subdivision,
const std::size_t z_subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Tet-Element mesh.

This algorithm generates regular grid points and split each grid cell into six tetrahedrals.

Parameters
x_lengthMesh dimension in x-direction.
y_lengthMesh dimension in y-direction.
z_lengthMesh dimension in z-direction.
x_subdivisionNumber of subdivisions in x-direction.
y_subdivisionNumber of subdivisions in y-direction.
z_subdivisionNumber of subdivisions in z-direction.
originOptional mesh's origin (the bottom lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 620 of file MeshGenerator.cpp.

629{
631 x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
632 y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
633}
MeshLib::Mesh * generateRegularTetMesh(const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, const BaseLib::ISubdivision &div_z, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")

References generateRegularTetMesh().

◆ generateRegularTetMesh() [3/3]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTetMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const unsigned n_z_cells,
const double cell_size_x,
const double cell_size_y,
const double cell_size_z,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 3D Tet-Element mesh.

This algorithm generates regular grid points and split each grid cell into six tetrahedrals.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
n_z_cellsNumber of cells in z-direction.
cell_size_xEdge length of Tet elements in x-direction.
cell_size_yEdge length of Tet elements in y s-direction.
cell_size_zEdge length of Tet elements in z-direction
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 635 of file MeshGenerator.cpp.

644{
646 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
647 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells),
648 BaseLib::UniformSubdivision(n_z_cells * cell_size_z, n_z_cells), origin,
649 mesh_name);
650}

References generateRegularTetMesh().

◆ generateRegularTriMesh() [1/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTriMesh ( const BaseLib::ISubdivision & div_x,
const BaseLib::ISubdivision & div_y,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Triangle-Element mesh. The mesh is generated in the x-y-plane.

Parameters
div_xSubdivision operator in x direction
div_ySubdivision operator in y direction
originOptional mesh's origin (the bottom lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 530 of file MeshGenerator.cpp.

535{
536 std::vector<double> vec_x(div_x());
537 std::vector<double> vec_y(div_y());
538 std::vector<MeshLib::Node*> nodes(
539 generateRegularNodes(vec_x, vec_y, origin));
540 const unsigned n_x_nodes(vec_x.size());
541 const unsigned n_x_cells(vec_x.size() - 1);
542 const unsigned n_y_cells(vec_y.size() - 1);
543
544 // elements
545 std::vector<MeshLib::Element*> elements;
546 elements.reserve(n_x_cells * n_y_cells * 2);
547
548 for (std::size_t j = 0; j < n_y_cells; j++)
549 {
550 const std::size_t offset_y1 = j * n_x_nodes;
551 const std::size_t offset_y2 = (j + 1) * n_x_nodes;
552 for (std::size_t k = 0; k < n_x_cells; k++)
553 {
554 elements.push_back(new MeshLib::Tri({nodes[offset_y1 + k],
555 nodes[offset_y2 + k + 1],
556 nodes[offset_y2 + k]}));
557
558 elements.push_back(new MeshLib::Tri({nodes[offset_y1 + k],
559 nodes[offset_y1 + k + 1],
560 nodes[offset_y2 + k + 1]}));
561 }
562 }
563
564 return new MeshLib::Mesh(mesh_name, nodes, elements,
565 true /* compute_element_neighbors */);
566}

References generateRegularNodes().

Referenced by CreateStructuredGridDialog::accept(), MeshToolsLib::RasterToMesh::convert(), generateRegularPrismMesh(), generateRegularTriMesh(), generateRegularTriMesh(), generateRegularTriMesh(), generateRegularTriMesh(), and main().

◆ generateRegularTriMesh() [2/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTriMesh ( const double length,
const std::size_t subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Triangle-Element mesh. The mesh is generated in the x-y-plane.

Parameters
lengthMesh's dimensions in x- and y-directions.
subdivisionNumber of subdivisions.
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 482 of file MeshGenerator.cpp.

487{
488 return generateRegularTriMesh(subdivision, subdivision,
489 length / subdivision, origin, mesh_name);
490}

References generateRegularTriMesh().

◆ generateRegularTriMesh() [3/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTriMesh ( const double x_length,
const double y_length,
const std::size_t x_subdivision,
const std::size_t y_subdivision,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Triangle-Element mesh. The mesh is generated in the x-y-plane.

Parameters
x_lengthMesh's dimension in x-direction.
y_lengthMesh's dimension in y-direction.
x_subdivisionNumber of subdivisions in x-direction.
y_subdivisionNumber of subdivisions in y-direction.
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 492 of file MeshGenerator.cpp.

499{
500 return generateRegularTriMesh(x_subdivision, y_subdivision,
501 x_length / x_subdivision,
502 y_length / y_subdivision, origin, mesh_name);
503}

References generateRegularTriMesh().

◆ generateRegularTriMesh() [4/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTriMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const double cell_size,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Triangle-Element mesh. The mesh is generated in the x-y-plane.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
cell_sizeEdge length of two equal sides of isosceles triangles
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 505 of file MeshGenerator.cpp.

511{
512 return generateRegularTriMesh(n_x_cells, n_y_cells, cell_size, cell_size,
513 origin, mesh_name);
514}

References generateRegularTriMesh().

◆ generateRegularTriMesh() [5/5]

MeshLib::Mesh * MeshToolsLib::MeshGenerator::generateRegularTriMesh ( const unsigned n_x_cells,
const unsigned n_y_cells,
const double cell_size_x,
const double cell_size_y,
MathLib::Point3d const & origin = MathLib::ORIGIN,
std::string const & mesh_name = "mesh" )

Generate a regular 2D Triangle-Element mesh. The mesh is generated in the x-y-plane.

Parameters
n_x_cellsNumber of cells in x-direction.
n_y_cellsNumber of cells in y-direction.
cell_size_xEdge length of triangles in x-direction.
cell_size_yEdge length of triangles in y-direction.
originOptional mesh's origin (the lower left corner) with MathLib::ORIGIN default.
mesh_nameName of the new mesh.

Definition at line 516 of file MeshGenerator.cpp.

523{
525 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
526 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells), origin,
527 mesh_name);
528}

References generateRegularTriMesh().