13 #include "libmesh/replicated_mesh.h"
14 #include "libmesh/mesh_generation.h"
15 #include "libmesh/string_to_enum.h"
16 #include "libmesh/periodic_boundaries.h"
17 #include "libmesh/periodic_boundary_base.h"
18 #include "libmesh/unstructured_mesh.h"
19 #include "libmesh/elem.h"
34 "EDGE EDGE2 EDGE3 EDGE4 QUAD QUAD4 QUAD8 QUAD9 TRI3 TRI6 HEX HEX8 HEX20 HEX27 TET4 TET10 "
35 "PRISM6 PRISM15 PRISM18 PYRAMID5 PYRAMID13 PYRAMID14");
40 params.
addParam<
unsigned int>(
"nx", 1,
"Number of elements in the X direction");
41 params.
addParam<
unsigned int>(
"ny", 1,
"Number of elements in the Y direction");
42 params.
addParam<
unsigned int>(
"nz", 1,
"Number of elements in the Z direction");
43 params.
addParam<Real>(
"xmin", 0.0,
"Lower X Coordinate of the generated mesh");
44 params.
addParam<Real>(
"ymin", 0.0,
"Lower Y Coordinate of the generated mesh");
45 params.
addParam<Real>(
"zmin", 0.0,
"Lower Z Coordinate of the generated mesh");
46 params.
addParam<Real>(
"xmax", 1.0,
"Upper X Coordinate of the generated mesh");
47 params.
addParam<Real>(
"ymax", 1.0,
"Upper Y Coordinate of the generated mesh");
48 params.
addParam<Real>(
"zmax", 1.0,
"Upper Z Coordinate of the generated mesh");
51 "The type of element from libMesh to "
52 "generate (default: linear element for "
53 "requested dimension)");
57 "Grade mesh into boundaries according to Gauss-Lobatto quadrature spacing.");
61 "bias_x>=0.5 & bias_x<=2",
62 "The amount by which to grow (or shrink) the cells in the x-direction.");
66 "bias_y>=0.5 & bias_y<=2",
67 "The amount by which to grow (or shrink) the cells in the y-direction.");
71 "bias_z>=0.5 & bias_z<=2",
72 "The amount by which to grow (or shrink) the cells in the z-direction.");
76 params.
addParam<std::vector<ExtraElementIDName>>(
"extra_element_integers",
77 "Names of extra element integers");
80 "Create a line, square, or cube mesh with uniformly spaced or biased elements.");
88 _nx(declareMeshProperty(
"num_elements_x", getParam<unsigned int>(
"nx"))),
89 _ny(declareMeshProperty(
"num_elements_y", getParam<unsigned int>(
"ny"))),
90 _nz(declareMeshProperty(
"num_elements_z", getParam<unsigned int>(
"nz"))),
91 _xmin(declareMeshProperty(
"xmin", getParam<Real>(
"xmin"))),
92 _xmax(declareMeshProperty(
"xmax", getParam<Real>(
"xmax"))),
93 _ymin(declareMeshProperty(
"ymin", getParam<Real>(
"ymin"))),
94 _ymax(declareMeshProperty(
"ymax", getParam<Real>(
"ymax"))),
95 _zmin(declareMeshProperty(
"zmin", getParam<Real>(
"zmin"))),
96 _zmax(declareMeshProperty(
"zmax", getParam<Real>(
"zmax"))),
97 _gauss_lobatto_grid(getParam<bool>(
"gauss_lobatto_grid")),
98 _bias_x(getParam<Real>(
"bias_x")),
99 _bias_y(getParam<Real>(
"bias_y")),
100 _bias_z(getParam<Real>(
"bias_z"))
103 mooseError(
"Cannot apply both Gauss-Lobatto mesh grading and biasing at the same time.");
106 std::unique_ptr<MeshBase>
110 auto mesh =
_mesh->buildMeshBaseObject();
114 for (
auto &
name : getParam<std::vector<ExtraElementIDName>>(
"extra_element_integers"))
115 mesh->add_elem_integer(
name);
118 MooseEnum elem_type_enum = getParam<MooseEnum>(
"elem_type");
126 elem_type_enum =
"EDGE2";
129 elem_type_enum =
"QUAD4";
132 elem_type_enum =
"HEX8";
137 ElemType elem_type = Utility::string_to_enum<ElemType>(elem_type_enum);
145 MeshTools::Generation::build_line(static_cast<UnstructuredMesh &>(*mesh),
153 MeshTools::Generation::build_square(static_cast<UnstructuredMesh &>(*mesh),
182 const auto MIN = std::numeric_limits<Real>::max();
185 std::array<Real, LIBMESH_DIM> bias = {
189 std::array<Real, LIBMESH_DIM> width = {
196 std::array<unsigned int, LIBMESH_DIM> nelem = {{
_nx,
_dim > 1 ?
_ny : 1,
_dim > 2 ?
_nz : 1}};
200 std::array<std::vector<Real>, LIBMESH_DIM> pows;
201 for (
unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
203 pows[dir].resize(nelem[dir] + 1);
204 for (
unsigned int i = 0; i < pows[dir].size(); ++i)
205 pows[dir][i] =
std::pow(bias[dir], static_cast<int>(i));
209 for (
auto & node_ptr : mesh->node_ptr_range())
211 Node & node = *node_ptr;
213 for (
unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
215 if (width[dir] != 0. && bias[dir] != 1.)
220 Real float_index = (node(dir) - mins[dir]) * nelem[dir] / width[dir];
222 Real integer_part = 0;
223 Real fractional_part = std::modf(float_index, &integer_part);
226 if (
std::abs(fractional_part) < TOLERANCE ||
std::abs(fractional_part - 1.0) < TOLERANCE)
235 int index =
round(float_index);
237 mooseAssert(index >= static_cast<int>(0) && index < static_cast<int>(pows[dir].size()),
238 "Scaled \"index\" out of range");
242 mins[dir] + width[dir] * (1. - pows[dir][index]) / (1. - pows[dir][nelem[dir]]);
244 else if (
std::abs(fractional_part - 0.5) < TOLERANCE)
255 node(dir) = mins[dir] +
257 (1. - 0.5 * (pows[dir][integer_part] + pows[dir][integer_part + 1])) /
258 (1. - pows[dir][nelem[dir]]);
263 mooseError(
"Unable to bias node at node(", dir,
")=", node(dir));
270 return dynamic_pointer_cast<MeshBase>(mesh);