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/node.h" 37 params.
addParam<
unsigned int>(
"nx", 1,
"Number of elements in the X direction");
38 params.
addParam<
unsigned int>(
"ny", 1,
"Number of elements in the Y direction");
39 params.
addParam<
unsigned int>(
"nz", 1,
"Number of elements in the Z direction");
40 params.
addParam<
Real>(
"xmin", 0.0,
"Lower X Coordinate of the generated mesh");
41 params.
addParam<
Real>(
"ymin", 0.0,
"Lower Y Coordinate of the generated mesh");
42 params.
addParam<
Real>(
"zmin", 0.0,
"Lower Z Coordinate of the generated mesh");
43 params.
addParam<
Real>(
"xmax", 1.0,
"Upper X Coordinate of the generated mesh");
44 params.
addParam<
Real>(
"ymax", 1.0,
"Upper Y Coordinate of the generated mesh");
45 params.
addParam<
Real>(
"zmax", 1.0,
"Upper Z Coordinate of the generated mesh");
48 "The type of element from libMesh to " 49 "generate (default: linear element for " 50 "requested dimension)");
54 "Grade mesh into boundaries according to Gauss-Lobatto quadrature spacing.");
58 "bias_x>=0.5 & bias_x<=2",
59 "The amount by which to grow (or shrink) the cells in the x-direction.");
63 "bias_y>=0.5 & bias_y<=2",
64 "The amount by which to grow (or shrink) the cells in the y-direction.");
68 "bias_z>=0.5 & bias_z<=2",
69 "The amount by which to grow (or shrink) the cells in the z-direction.");
74 "Create a line, square, or cube mesh with uniformly spaced or biased elements.");
81 _nx(getParam<unsigned
int>(
"nx")),
82 _ny(getParam<unsigned
int>(
"ny")),
83 _nz(getParam<unsigned
int>(
"nz")),
84 _xmin(getParam<
Real>(
"xmin")),
85 _xmax(getParam<
Real>(
"xmax")),
86 _ymin(getParam<
Real>(
"ymin")),
87 _ymax(getParam<
Real>(
"ymax")),
88 _zmin(getParam<
Real>(
"zmin")),
89 _zmax(getParam<
Real>(
"zmax")),
90 _gauss_lobatto_grid(getParam<bool>(
"gauss_lobatto_grid")),
91 _bias_x(getParam<
Real>(
"bias_x")),
92 _bias_y(getParam<
Real>(
"bias_y")),
93 _bias_z(getParam<
Real>(
"bias_z")),
94 _dims_may_have_changed(false)
97 mooseError(
"Cannot apply both Gauss-Lobatto mesh grading and biasing at the same time.");
151 std::unique_ptr<MooseMesh>
160 MooseEnum elem_type_enum = getParam<MooseEnum>(
"elem_type");
168 elem_type_enum =
"EDGE2";
171 elem_type_enum =
"QUAD4";
174 elem_type_enum =
"HEX8";
179 ElemType elem_type = Utility::string_to_enum<ElemType>(elem_type_enum);
187 MeshTools::Generation::build_line(dynamic_cast<UnstructuredMesh &>(
getMesh()),
195 MeshTools::Generation::build_square(dynamic_cast<UnstructuredMesh &>(
getMesh()),
206 MeshTools::Generation::build_cube(dynamic_cast<UnstructuredMesh &>(
getMesh()),
228 std::array<Real, LIBMESH_DIM> bias = {
232 std::array<Real, LIBMESH_DIM> width = {
236 std::array<Real, LIBMESH_DIM> mins = {
240 std::array<unsigned int, LIBMESH_DIM> nelem = {{
_nx,
_dim > 1 ?
_ny : 1,
_dim > 2 ?
_nz : 1}};
244 std::array<std::vector<Real>, LIBMESH_DIM> pows;
245 for (
unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
247 pows[dir].resize(nelem[dir] + 1);
248 for (
unsigned int i = 0; i < pows[dir].size(); ++i)
249 pows[dir][i] =
std::pow(bias[dir], static_cast<int>(i));
253 for (
auto & node_ptr :
mesh.node_ptr_range())
255 Node &
node = *node_ptr;
257 for (
unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
259 if (width[dir] != 0. && bias[dir] != 1.)
264 Real float_index = (
node(dir) - mins[dir]) * nelem[dir] / width[dir];
266 Real integer_part = 0;
267 Real fractional_part = std::modf(float_index, &integer_part);
270 if (
std::abs(fractional_part) < TOLERANCE ||
std::abs(fractional_part - 1.0) < TOLERANCE)
279 int index =
round(float_index);
281 mooseAssert(index >= static_cast<int>(0) && index < static_cast<int>(pows[dir].size()),
282 "Scaled \"index\" out of range");
286 mins[dir] + width[dir] * (1. - pows[dir][index]) / (1. - pows[dir][nelem[dir]]);
288 else if (
std::abs(fractional_part - 0.5) < TOLERANCE)
299 node(dir) = mins[dir] +
301 (1. - 0.5 * (pows[dir][integer_part] + pows[dir][integer_part + 1])) /
302 (1. - pows[dir][nelem[dir]]);
static InputParameters validParams()
Typical "Moose-style" constructor and copy constructor.
virtual Real getMaxInDimension(unsigned int component) const
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
virtual Real getMinInDimension(unsigned int component) const
Returns the min or max of the requested dimension respectively.
bool prepared() const
Setter/getter for whether the mesh is prepared.
const std::string LIST_GEOM_ELEM
virtual std::unique_ptr< MooseMesh > safeClone() const override
A safer version of the clone() method that hands back an allocated object wrapped in a smart pointer...
MooseEnum _dim
The dimension of the mesh.
virtual Real getMaxInDimension(unsigned int component) const override
bool _dims_may_have_changed
Boolean to indicate that dimensions may have changed.
virtual Real getMinInDimension(unsigned int component) const override
Returns the min or max of the requested dimension respectively.
std::unique_ptr< T > copyConstruct(const T &object)
Copy constructs the object object.
virtual void buildMesh() override
Must be overridden by child classes.
Factory & getFactory()
Retrieve a writable reference to the Factory associated with this App.
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Mesh generated from parameters.
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Real _xmin
The min/max values for x,y,z component.
GeneratedMesh(const InputParameters ¶meters)
MooseApp & _app
The MOOSE application this is associated with.
registerMooseObject("MooseApp", GeneratedMesh)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _gauss_lobatto_grid
All of the libmesh build_line/square/cube routines support an option to grade the mesh into the bound...
virtual const Node & node(const dof_id_type i) const
Various accessors (pointers/references) for Node "i".
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
unsigned int _nx
Number of elements in x, y, z direction.
Real _bias_x
The amount by which to bias the cells in the x,y,z directions.
MooseUnits pow(const MooseUnits &, int)
void ErrorVector unsigned int
bool _regular_orthogonal_mesh
Boolean indicating whether this mesh was detected to be regular and orthogonal.
static InputParameters validParams()