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 747 of file MeshGenerator.cpp.

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

◆ 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 179 of file MeshGenerator.cpp.

182{
183 const std::vector<double> vec_x(div());
184 std::vector<MeshLib::Node*> nodes(generateRegularNodes(vec_x, origin));
185
186 // elements
187 const std::size_t n_cells = nodes.size() - 1;
188 std::vector<MeshLib::Element*> elements;
189 elements.reserve(n_cells);
190
191 for (std::size_t i = 0; i < n_cells; i++)
192 {
193 elements.push_back(new MeshLib::Line({nodes[i], nodes[i + 1]}));
194 }
195
196 return new MeshLib::Mesh(mesh_name, nodes, elements,
197 true /* compute_element_neighbors */);
198}
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 160 of file MeshGenerator.cpp.

164{
165 return generateLineMesh(subdivision, length / subdivision, origin,
166 mesh_name);
167}
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 169 of file MeshGenerator.cpp.

173{
174 return generateLineMesh(
175 BaseLib::UniformSubdivision(n_cells * cell_size, n_cells), origin,
176 mesh_name);
177}

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 339 of file MeshGenerator.cpp.

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

References generateRegularNodes().

Referenced by CreateStructuredGridDialog::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 283 of file MeshGenerator.cpp.

288{
289 return MeshGenerator::generateRegularHexMesh(
290 subdivision, subdivision, subdivision, length / subdivision, origin,
291 mesh_name);
292}

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 294 of file MeshGenerator.cpp.

303{
304 return MeshGenerator::generateRegularHexMesh(
305 x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
306 y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
307}

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 309 of file MeshGenerator.cpp.

316{
317 return MeshGenerator::generateRegularHexMesh(
318 n_x_cells, n_y_cells, n_z_cells, cell_size, cell_size, cell_size,
319 origin, mesh_name);
320}

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 322 of file MeshGenerator.cpp.

331{
333 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
334 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells),
335 BaseLib::UniformSubdivision(n_z_cells * cell_size_z, n_z_cells), origin,
336 mesh_name);
337}
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().

◆ 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 103 of file MeshGenerator.cpp.

107{
108 std::vector<MeshLib::Node*> nodes;
109 nodes.reserve((n_cells[0] + 1) * (n_cells[1] + 1) * (n_cells[2] + 1));
110
111 for (std::size_t i = 0; i < n_cells[2] + 1; i++)
112 {
113 const double z(origin[2] + cell_size[2] * i);
114 for (std::size_t j = 0; j < n_cells[1] + 1; j++)
115 {
116 const double y(origin[1] + cell_size[1] * j);
117 for (std::size_t k = 0; k < n_cells[0] + 1; k++)
118 {
119 nodes.push_back(
120 new MeshLib::Node(origin[0] + cell_size[0] * k, y, z));
121 }
122 }
123 }
124 return nodes;
125}

◆ 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 28 of file MeshGenerator.cpp.

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

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 61 of file MeshGenerator.cpp.

63{
64 std::vector<const std::vector<double>*> vec_xyz_coords;
65 vec_xyz_coords.push_back(&vec_x_coords);
66 std::vector<double> dummy(1, 0.0);
67 for (unsigned i = vec_xyz_coords.size() - 1; i < 3u; i++)
68 {
69 vec_xyz_coords.push_back(&dummy);
70 }
71 return generateRegularNodes(vec_xyz_coords, origin);
72}

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 74 of file MeshGenerator.cpp.

78{
79 std::vector<const std::vector<double>*> vec_xyz_coords;
80 vec_xyz_coords.push_back(&vec_x_coords);
81 vec_xyz_coords.push_back(&vec_y_coords);
82 std::vector<double> dummy(1, 0.0);
83 for (unsigned i = vec_xyz_coords.size() - 1; i < 3u; i++)
84 {
85 vec_xyz_coords.push_back(&dummy);
86 }
87 return generateRegularNodes(vec_xyz_coords, origin);
88}

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 90 of file MeshGenerator.cpp.

95{
96 std::vector<const std::vector<double>*> vec_xyz_coords;
97 vec_xyz_coords.push_back(&vec_x_coords);
98 vec_xyz_coords.push_back(&vec_y_coords);
99 vec_xyz_coords.push_back(&vec_z_coords);
100 return generateRegularNodes(vec_xyz_coords, origin);
101}

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 575 of file MeshGenerator.cpp.

584{
586 x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
587 y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
588}
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 590 of file MeshGenerator.cpp.

597{
598 return generateRegularPrismMesh(n_x_cells, n_y_cells, n_z_cells, cell_size,
599 cell_size, cell_size, origin, mesh_name);
600}

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 602 of file MeshGenerator.cpp.

611{
612 std::unique_ptr<MeshLib::Mesh> mesh(generateRegularTriMesh(
613 n_x_cells, n_y_cells, cell_size_x, cell_size_y, origin, mesh_name));
614 std::size_t const n_tris(mesh->getNumberOfElements());
615 bool const copy_material_ids = false;
616 for (std::size_t i = 0; i < n_z_cells; ++i)
617 {
618 mesh.reset(MeshToolsLib::addLayerToMesh(*mesh, cell_size_z, mesh_name,
619 true, copy_material_ids));
620 }
621 std::vector<std::size_t> elem_ids(n_tris);
622 std::iota(elem_ids.begin(), elem_ids.end(), 0);
623 return MeshToolsLib::removeElements(*mesh, elem_ids, mesh_name);
624}
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 thickness, std::string const &name, bool on_top, bool copy_material_ids)
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 391 of file MeshGenerator.cpp.

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

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

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 249 of file MeshGenerator.cpp.

254{
255 std::vector<double> vec_x(div_x());
256 std::vector<double> vec_y(div_y());
257 std::vector<MeshLib::Node*> nodes(
258 generateRegularNodes(vec_x, vec_y, origin));
259 const unsigned n_x_nodes(vec_x.size());
260
261 // elements
262 std::vector<MeshLib::Element*> elements;
263 const unsigned n_x_cells(vec_x.size() - 1);
264 const unsigned n_y_cells(vec_y.size() - 1);
265 elements.reserve(n_x_cells * n_y_cells);
266
267 for (std::size_t j = 0; j < n_y_cells; j++)
268 {
269 const std::size_t offset_y1 = j * n_x_nodes;
270 const std::size_t offset_y2 = (j + 1) * n_x_nodes;
271 for (std::size_t k = 0; k < n_x_cells; k++)
272 {
273 elements.push_back(new MeshLib::Quad(
274 {nodes[offset_y1 + k], nodes[offset_y1 + k + 1],
275 nodes[offset_y2 + k + 1], nodes[offset_y2 + k]}));
276 }
277 }
278
279 return new MeshLib::Mesh(mesh_name, nodes, elements,
280 true /* compute_element_neighbors */);
281}

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 200 of file MeshGenerator.cpp.

205{
206 return generateRegularQuadMesh(subdivision, subdivision,
207 length / subdivision, length / subdivision,
208 origin, mesh_name);
209}
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 211 of file MeshGenerator.cpp.

218{
219 return generateRegularQuadMesh(x_subdivision, y_subdivision,
220 x_length / x_subdivision,
221 y_length / y_subdivision, origin, mesh_name);
222}

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 224 of file MeshGenerator.cpp.

230{
231 return generateRegularQuadMesh(n_x_cells, n_y_cells, cell_size, cell_size,
232 origin, mesh_name);
233}

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 235 of file MeshGenerator.cpp.

242{
244 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
245 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells), origin,
246 mesh_name);
247}

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 658 of file MeshGenerator.cpp.

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

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 626 of file MeshGenerator.cpp.

635{
637 x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
638 y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
639}
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 641 of file MeshGenerator.cpp.

650{
652 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
653 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells),
654 BaseLib::UniformSubdivision(n_z_cells * cell_size_z, n_z_cells), origin,
655 mesh_name);
656}

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 537 of file MeshGenerator.cpp.

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

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 489 of file MeshGenerator.cpp.

494{
495 return generateRegularTriMesh(subdivision, subdivision,
496 length / subdivision, origin, mesh_name);
497}

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 499 of file MeshGenerator.cpp.

506{
507 return generateRegularTriMesh(x_subdivision, y_subdivision,
508 x_length / x_subdivision,
509 y_length / y_subdivision, origin, mesh_name);
510}

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 512 of file MeshGenerator.cpp.

518{
519 return generateRegularTriMesh(n_x_cells, n_y_cells, cell_size, cell_size,
520 origin, mesh_name);
521}

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 523 of file MeshGenerator.cpp.

530{
532 BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
533 BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells), origin,
534 mesh_name);
535}

References generateRegularTriMesh().