OGS
MeshLib::MeshGenerator Namespace Reference

Detailed Description

Functionality for generating meshing.

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)
 
MeshgenerateLineMesh (const BaseLib::ISubdivision &div, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
MeshgenerateLineMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
MeshgenerateLineMesh (const unsigned n_cells, const double cell_size, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
MeshgenerateRegularQuadMesh (const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
MeshgenerateRegularQuadMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
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")
 
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")
 
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")
 
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")
 
MeshgenerateRegularHexMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
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")
 
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")
 
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")
 
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")
 
MeshgenerateRegularTriMesh (const BaseLib::ISubdivision &div_x, const BaseLib::ISubdivision &div_y, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
MeshgenerateRegularTriMesh (const double length, const std::size_t subdivision, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")
 
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")
 
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")
 
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")
 
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")
 
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")
 
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")
 
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")
 
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")
 
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 * MeshLib::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 713 of file MeshGenerator.cpp.

717 {
718  std::array<double, 2> step_size{{(ur[0] - ll[0]) / (n_steps[0] - 1),
719  (ur[1] - ll[1]) / (n_steps[1] - 1)}};
720 
721  std::vector<MeshLib::Node*> nodes;
722  for (std::size_t j(0); j < n_steps[1]; ++j)
723  {
724  for (std::size_t i(0); i < n_steps[0]; ++i)
725  {
726  std::size_t const id = i + j * n_steps[1];
727  double const x = ll[0] + i * step_size[0];
728  double const y = ll[1] + j * step_size[1];
729 
730  nodes.push_back(new MeshLib::Node({x, y, f(x, y)}, id));
731  }
732  }
733 
734  std::vector<MeshLib::Element*> sfc_eles;
735  for (std::size_t j(0); j < n_steps[1] - 1; ++j)
736  {
737  for (std::size_t i(0); i < n_steps[0] - 1; ++i)
738  {
739  std::size_t id_ll(i + j * n_steps[0]);
740  std::size_t id_lr(i + 1 + j * n_steps[0]);
741  std::size_t id_ul(i + (j + 1) * n_steps[0]);
742  std::size_t id_ur(i + 1 + (j + 1) * n_steps[0]);
743  sfc_eles.push_back(
744  new MeshLib::Tri({nodes[id_ll], nodes[id_lr], nodes[id_ur]}));
745  sfc_eles.push_back(
746  new MeshLib::Tri({nodes[id_ll], nodes[id_ur], nodes[id_ul]}));
747  }
748  }
749 
750  return new MeshLib::Mesh(mesh_name, nodes, sfc_eles);
751 }

◆ generateLineMesh() [1/3]

Mesh * MeshLib::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 177 of file MeshGenerator.cpp.

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

References generateRegularNodes().

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

◆ generateLineMesh() [2/3]

Mesh * MeshLib::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 158 of file MeshGenerator.cpp.

162 {
163  return generateLineMesh(subdivision, length / subdivision, origin,
164  mesh_name);
165 }
Mesh * generateLineMesh(const BaseLib::ISubdivision &div, MathLib::Point3d const &origin=MathLib::ORIGIN, std::string const &mesh_name="mesh")

References generateLineMesh().

◆ generateLineMesh() [3/3]

Mesh * MeshLib::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 167 of file MeshGenerator.cpp.

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

References generateLineMesh().

◆ generateRegularHexMesh() [1/5]

Mesh * MeshLib::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 325 of file MeshGenerator.cpp.

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

References generateRegularNodes().

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

◆ generateRegularHexMesh() [2/5]

Mesh * MeshLib::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 273 of file MeshGenerator.cpp.

277 {
279  subdivision, subdivision, subdivision, length / subdivision, origin,
280  mesh_name);
281 }
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]

Mesh * MeshLib::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 283 of file MeshGenerator.cpp.

291 {
293  x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
294  y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
295 }

References generateRegularHexMesh().

◆ generateRegularHexMesh() [4/5]

Mesh * MeshLib::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 297 of file MeshGenerator.cpp.

303 {
305  n_x_cells, n_y_cells, n_z_cells, cell_size, cell_size, cell_size,
306  origin, mesh_name);
307 }

References generateRegularHexMesh().

◆ generateRegularHexMesh() [5/5]

Mesh * MeshLib::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 309 of file MeshGenerator.cpp.

317 {
318  return generateRegularHexMesh(
319  BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
320  BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells),
321  BaseLib::UniformSubdivision(n_z_cells * cell_size_z, n_z_cells), origin,
322  mesh_name);
323 }

References generateRegularHexMesh().

◆ generateRegularNodes() [1/5]

std::vector< MeshLib::Node * > MeshLib::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 102 of file MeshGenerator.cpp.

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

◆ generateRegularNodes() [2/5]

std::vector< MeshLib::Node * > MeshLib::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<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(
53  coords[0].begin(), coords[0].end(), std::back_inserter(nodes),
54  [&y, &z](double const& x) { return new Node(x, y, z); });
55  }
56  }
57  return nodes;
58 }

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

◆ generateRegularNodes() [3/5]

std::vector< MeshLib::Node * > MeshLib::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 60 of file MeshGenerator.cpp.

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

References generateRegularNodes().

◆ generateRegularNodes() [4/5]

std::vector< MeshLib::Node * > MeshLib::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 73 of file MeshGenerator.cpp.

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

References generateRegularNodes().

◆ generateRegularNodes() [5/5]

std::vector< MeshLib::Node * > MeshLib::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 89 of file MeshGenerator.cpp.

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

References generateRegularNodes().

◆ generateRegularPrismMesh() [1/3]

Mesh * MeshLib::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 549 of file MeshGenerator.cpp.

557 {
559  x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
560  y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
561 }
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")

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

◆ generateRegularPrismMesh() [2/3]

Mesh * MeshLib::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 563 of file MeshGenerator.cpp.

569 {
570  return generateRegularPrismMesh(n_x_cells, n_y_cells, n_z_cells, cell_size,
571  cell_size, cell_size, origin, mesh_name);
572 }

References generateRegularPrismMesh().

◆ generateRegularPrismMesh() [3/3]

Mesh * MeshLib::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 574 of file MeshGenerator.cpp.

582 {
583  std::unique_ptr<MeshLib::Mesh> mesh(generateRegularTriMesh(
584  n_x_cells, n_y_cells, cell_size_x, cell_size_y, origin, mesh_name));
585  std::size_t const n_tris(mesh->getNumberOfElements());
586  bool const copy_material_ids = false;
587  for (std::size_t i = 0; i < n_z_cells; ++i)
588  {
589  mesh.reset(MeshLib::addLayerToMesh(*mesh, cell_size_z, mesh_name, true,
590  copy_material_ids));
591  }
592  std::vector<std::size_t> elem_ids(n_tris);
593  std::iota(elem_ids.begin(), elem_ids.end(), 0);
594  return MeshLib::removeElements(*mesh, elem_ids, mesh_name);
595 }
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 * removeElements(const MeshLib::Mesh &mesh, const std::vector< std::size_t > &removed_element_ids, const std::string &new_mesh_name)
MeshLib::Mesh * addLayerToMesh(MeshLib::Mesh const &mesh, double thickness, std::string const &name, bool on_top, bool copy_material_ids)

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

◆ generateRegularPyramidMesh()

Mesh * MeshLib::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 374 of file MeshGenerator.cpp.

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

References generateRegularNodes(), and generateRegularPyramidTopNodes().

Referenced by main().

◆ generateRegularPyramidTopNodes()

std::vector<MeshLib::Node*> MeshLib::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 127 of file MeshGenerator.cpp.

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

Referenced by generateRegularPyramidMesh().

◆ generateRegularQuadMesh() [1/5]

Mesh * MeshLib::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.

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

References generateRegularNodes().

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

◆ generateRegularQuadMesh() [2/5]

Mesh * MeshLib::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 197 of file MeshGenerator.cpp.

201 {
202  return generateRegularQuadMesh(subdivision, subdivision,
203  length / subdivision, length / subdivision,
204  origin, mesh_name);
205 }
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]

Mesh * MeshLib::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 207 of file MeshGenerator.cpp.

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

References generateRegularQuadMesh().

◆ generateRegularQuadMesh() [4/5]

Mesh * MeshLib::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 219 of file MeshGenerator.cpp.

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

References generateRegularQuadMesh().

◆ generateRegularQuadMesh() [5/5]

Mesh * MeshLib::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 229 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]

Mesh * MeshLib::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 627 of file MeshGenerator.cpp.

632 {
633  std::vector<double> vec_x(div_x());
634  std::vector<double> vec_y(div_y());
635  std::vector<double> vec_z(div_z());
636  std::vector<Node*> nodes(generateRegularNodes(vec_x, vec_y, vec_z, origin));
637 
638  const unsigned n_x_nodes(vec_x.size());
639  const unsigned n_y_nodes(vec_y.size());
640  const unsigned n_x_cells(vec_x.size() - 1);
641  const unsigned n_y_cells(vec_y.size() - 1);
642  const unsigned n_z_cells(vec_z.size() - 1);
643 
644  // elements
645  std::vector<Element*> elements;
646  elements.reserve(n_x_cells * n_y_cells * n_z_cells * 6);
647 
648  for (std::size_t i = 0; i < n_z_cells; i++)
649  {
650  const std::size_t offset_z1 = i * n_x_nodes * n_y_nodes; // bottom
651  const std::size_t offset_z2 = (i + 1) * n_x_nodes * n_y_nodes; // top
652  for (std::size_t j = 0; j < n_y_cells; j++)
653  {
654  const std::size_t offset_y1 = j * n_x_nodes;
655  const std::size_t offset_y2 = (j + 1) * n_x_nodes;
656  for (std::size_t k = 0; k < n_x_cells; k++)
657  {
658  // tet 1
659  elements.push_back(
660  new Tet({// bottom
661  nodes[offset_z1 + offset_y1 + k],
662  nodes[offset_z1 + offset_y2 + k + 1],
663  nodes[offset_z1 + offset_y2 + k],
664  // top
665  nodes[offset_z2 + offset_y1 + k]}));
666  // tet 2
667  elements.push_back(
668  new Tet({// bottom
669  nodes[offset_z1 + offset_y2 + k + 1],
670  nodes[offset_z1 + offset_y2 + k],
671  // top
672  nodes[offset_z2 + offset_y1 + k],
673  nodes[offset_z2 + offset_y2 + k + 1]}));
674  // tet 3
675  elements.push_back(
676  new Tet({// bottom
677  nodes[offset_z1 + offset_y2 + k],
678  // top
679  nodes[offset_z2 + offset_y1 + k],
680  nodes[offset_z2 + offset_y2 + k + 1],
681  nodes[offset_z2 + offset_y2 + k]}));
682  // tet 4
683  elements.push_back(
684  new Tet({// bottom
685  nodes[offset_z1 + offset_y1 + k],
686  nodes[offset_z1 + offset_y1 + k + 1],
687  nodes[offset_z1 + offset_y2 + k + 1],
688  // top
689  nodes[offset_z2 + offset_y1 + k + 1]}));
690  // tet 5
691  elements.push_back(
692  new Tet({// bottom
693  nodes[offset_z1 + offset_y1 + k],
694  nodes[offset_z1 + offset_y2 + k + 1],
695  // top
696  nodes[offset_z2 + offset_y1 + k],
697  nodes[offset_z2 + offset_y1 + k + 1]}));
698  // tet 6
699  elements.push_back(
700  new Tet({// bottom
701  nodes[offset_z1 + offset_y2 + k + 1],
702  // top
703  nodes[offset_z2 + offset_y1 + k],
704  nodes[offset_z2 + offset_y1 + k + 1],
705  nodes[offset_z2 + offset_y2 + k + 1]}));
706  }
707  }
708  }
709 
710  return new Mesh(mesh_name, nodes, elements);
711 }
TemplateElement< MeshLib::TetRule4 > Tet
Definition: Tet.h:25

References generateRegularNodes().

Referenced by generateRegularTetMesh(), and main().

◆ generateRegularTetMesh() [2/3]

Mesh * MeshLib::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 597 of file MeshGenerator.cpp.

605 {
606  return generateRegularTetMesh(
607  x_subdivision, y_subdivision, z_subdivision, x_length / x_subdivision,
608  y_length / y_subdivision, z_length / z_subdivision, origin, mesh_name);
609 }
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]

Mesh * MeshLib::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 611 of file MeshGenerator.cpp.

619 {
620  return generateRegularTetMesh(
621  BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
622  BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells),
623  BaseLib::UniformSubdivision(n_z_cells * cell_size_z, n_z_cells), origin,
624  mesh_name);
625 }

References generateRegularTetMesh().

◆ generateRegularTriMesh() [1/5]

Mesh * MeshLib::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 514 of file MeshGenerator.cpp.

518 {
519  std::vector<double> vec_x(div_x());
520  std::vector<double> vec_y(div_y());
521  std::vector<Node*> nodes(generateRegularNodes(vec_x, vec_y, origin));
522  const unsigned n_x_nodes(vec_x.size());
523  const unsigned n_x_cells(vec_x.size() - 1);
524  const unsigned n_y_cells(vec_y.size() - 1);
525 
526  // elements
527  std::vector<Element*> elements;
528  elements.reserve(n_x_cells * n_y_cells * 2);
529 
530  for (std::size_t j = 0; j < n_y_cells; j++)
531  {
532  const std::size_t offset_y1 = j * n_x_nodes;
533  const std::size_t offset_y2 = (j + 1) * n_x_nodes;
534  for (std::size_t k = 0; k < n_x_cells; k++)
535  {
536  elements.push_back(
537  new Tri({nodes[offset_y1 + k], nodes[offset_y2 + k + 1],
538  nodes[offset_y2 + k]}));
539 
540  elements.push_back(
541  new Tri({nodes[offset_y1 + k], nodes[offset_y1 + k + 1],
542  nodes[offset_y2 + k + 1]}));
543  }
544  }
545 
546  return new Mesh(mesh_name, nodes, elements);
547 }
TemplateElement< MeshLib::TriRule3 > Tri
Definition: Tri.h:26

References generateRegularNodes().

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

◆ generateRegularTriMesh() [2/5]

Mesh * MeshLib::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 470 of file MeshGenerator.cpp.

474 {
475  return generateRegularTriMesh(subdivision, subdivision,
476  length / subdivision, origin, mesh_name);
477 }

References generateRegularTriMesh().

◆ generateRegularTriMesh() [3/5]

Mesh * MeshLib::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 479 of file MeshGenerator.cpp.

485 {
486  return generateRegularTriMesh(x_subdivision, y_subdivision,
487  x_length / x_subdivision,
488  y_length / y_subdivision, origin, mesh_name);
489 }

References generateRegularTriMesh().

◆ generateRegularTriMesh() [4/5]

Mesh * MeshLib::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 491 of file MeshGenerator.cpp.

496 {
497  return generateRegularTriMesh(n_x_cells, n_y_cells, cell_size, cell_size,
498  origin, mesh_name);
499 }

References generateRegularTriMesh().

◆ generateRegularTriMesh() [5/5]

Mesh * MeshLib::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 501 of file MeshGenerator.cpp.

507 {
508  return generateRegularTriMesh(
509  BaseLib::UniformSubdivision(n_x_cells * cell_size_x, n_x_cells),
510  BaseLib::UniformSubdivision(n_y_cells * cell_size_y, n_y_cells), origin,
511  mesh_name);
512 }

References generateRegularTriMesh().