13 #include "libmesh/replicated_mesh.h" 14 #include "libmesh/mesh_generation.h" 15 #include "libmesh/mesh_modification.h" 16 #include "libmesh/string_to_enum.h" 17 #include "libmesh/periodic_boundaries.h" 18 #include "libmesh/periodic_boundary_base.h" 19 #include "libmesh/unstructured_mesh.h" 20 #include "libmesh/elem.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)");
51 params.
addParam<std::vector<SubdomainID>>(
53 "Subdomain IDs for each element, default to all zero. If a single number is specified, that " 54 "subdomain id is used for all element.");
55 params.
addParam<SubdomainName>(
"subdomain_name",
56 "If specified, single subdomain name for all elements");
61 "Grade mesh into boundaries according to Gauss-Lobatto quadrature spacing.");
65 "bias_x>=0.5 & bias_x<=2",
66 "The amount by which to grow (or shrink) the cells in the x-direction.");
70 "bias_y>=0.5 & bias_y<=2",
71 "The amount by which to grow (or shrink) the cells in the y-direction.");
75 "bias_z>=0.5 & bias_z<=2",
76 "The amount by which to grow (or shrink) the cells in the z-direction.");
78 params.
addParam<std::string>(
"boundary_name_prefix",
79 "If provided, prefix the built in boundary names with this string");
81 "boundary_id_offset", 0,
"This offset is added to the generated boundary IDs");
83 params.
addParam<std::vector<ExtraElementIDName>>(
"extra_element_integers",
84 "Names of extra element integers");
87 "Create a line, square, or cube mesh with uniformly spaced or biased elements.");
95 _nx(declareMeshProperty(
"num_elements_x", getParam<unsigned
int>(
"nx"))),
96 _ny(declareMeshProperty(
"num_elements_y", getParam<unsigned
int>(
"ny"))),
97 _nz(declareMeshProperty(
"num_elements_z", getParam<unsigned
int>(
"nz"))),
98 _xmin(declareMeshProperty(
"xmin", getParam<
Real>(
"xmin"))),
99 _xmax(declareMeshProperty(
"xmax", getParam<
Real>(
"xmax"))),
100 _ymin(declareMeshProperty(
"ymin", getParam<
Real>(
"ymin"))),
101 _ymax(declareMeshProperty(
"ymax", getParam<
Real>(
"ymax"))),
102 _zmin(declareMeshProperty(
"zmin", getParam<
Real>(
"zmin"))),
103 _zmax(declareMeshProperty(
"zmax", getParam<
Real>(
"zmax"))),
104 _has_subdomain_ids(isParamValid(
"subdomain_ids")),
105 _gauss_lobatto_grid(getParam<bool>(
"gauss_lobatto_grid")),
106 _bias_x(getParam<
Real>(
"bias_x")),
107 _bias_y(getParam<
Real>(
"bias_y")),
108 _bias_z(getParam<
Real>(
"bias_z")),
109 _boundary_name_prefix(isParamValid(
"boundary_name_prefix")
110 ? getParam<
std::string>(
"boundary_name_prefix") +
"_" 115 mooseError(
"Cannot apply both Gauss-Lobatto mesh grading and biasing at the same time.");
117 paramError(
"xmax",
"xmax must be larger than xmin.");
119 paramError(
"ymax",
"ymax must be larger than ymin.");
121 paramError(
"zmax",
"zmax must be larger than zmin.");
124 std::unique_ptr<MeshBase>
132 for (
auto &
name :
getParam<std::vector<ExtraElementIDName>>(
"extra_element_integers"))
136 MooseEnum elem_type_enum = getParam<MooseEnum>(
"elem_type");
144 elem_type_enum =
"EDGE2";
147 elem_type_enum =
"QUAD4";
150 elem_type_enum =
"HEX8";
155 ElemType elem_type = Utility::string_to_enum<ElemType>(elem_type_enum);
163 MeshTools::Generation::build_line(static_cast<UnstructuredMesh &>(*
mesh),
171 MeshTools::Generation::build_square(static_cast<UnstructuredMesh &>(*
mesh),
182 MeshTools::Generation::build_cube(static_cast<UnstructuredMesh &>(*
mesh),
199 auto & bids = getParam<std::vector<SubdomainID>>(
"subdomain_ids");
200 if (bids.size() !=
_nx *
_ny *
_nz && bids.size() != 1)
202 "Size must equal to the product of number of elements in all directions, or one.");
203 for (
auto & elem :
mesh->element_ptr_range())
205 const Point p = elem->vertex_average();
209 unsigned int i = iz *
_nx *
_ny + iy *
_nx + ix;
210 if (bids.size() == 1)
211 elem->subdomain_id() = bids[0];
213 elem->subdomain_id() = bids[i];
219 const auto & subdomain_name = getParam<SubdomainName>(
"subdomain_name");
222 const auto & bids = getParam<std::vector<SubdomainID>>(
"subdomain_ids");
226 "Specifying a subdomain_name is only supported for a single entry in subdomain_ids");
228 mesh->subdomain_name(bids[0]) = subdomain_name;
231 mesh->subdomain_name(0) = subdomain_name;
235 BoundaryInfo & boundary_info =
mesh->get_boundary_info();
238 const auto mesh_boundary_ids = boundary_info.get_global_boundary_ids();
239 for (
auto rit = mesh_boundary_ids.rbegin(); rit != mesh_boundary_ids.rend(); ++rit)
241 const std::string old_sideset_name = boundary_info.sideset_name(*rit);
242 const std::string old_nodeset_name = boundary_info.nodeset_name(*rit);
258 std::array<Real, LIBMESH_DIM> bias = {
262 std::array<Real, LIBMESH_DIM> width = {
269 std::array<unsigned int, LIBMESH_DIM> nelem = {{
_nx,
_dim > 1 ?
_ny : 1,
_dim > 2 ?
_nz : 1}};
273 std::array<std::vector<Real>, LIBMESH_DIM> pows;
274 for (
unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
276 pows[dir].resize(nelem[dir] + 1);
278 for (
unsigned int i = 1; i < pows[dir].size(); ++i)
279 pows[dir][i] = pows[dir][i - 1] * bias[dir];
283 for (
auto & node_ptr :
mesh->node_ptr_range())
285 Node & node = *node_ptr;
287 for (
unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
289 if (width[dir] != 0. && bias[dir] != 1.)
294 Real float_index = (node(dir) - mins[dir]) * nelem[dir] / width[dir];
296 Real integer_part = 0;
297 Real fractional_part = std::modf(float_index, &integer_part);
300 if (
std::abs(fractional_part) < TOLERANCE ||
std::abs(fractional_part - 1.0) < TOLERANCE)
309 int index =
round(float_index);
311 mooseAssert(index >= static_cast<int>(0) && index < static_cast<int>(pows[dir].size()),
312 "Scaled \"index\" out of range");
316 mins[dir] + width[dir] * (1. - pows[dir][index]) / (1. - pows[dir][nelem[dir]]);
318 else if (
std::abs(fractional_part - 0.5) < TOLERANCE)
329 node(dir) = mins[dir] +
331 (1. - 0.5 * (pows[dir][integer_part] + pows[dir][integer_part + 1])) /
332 (1. - pows[dir][nelem[dir]]);
337 mooseError(
"Unable to bias node at node(", dir,
")=", node(dir));
344 mesh->set_isnt_prepared();
const Real _bias_x
The amount by which to bias the cells in the x,y,z directions.
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
const std::string LIST_GEOM_ELEM
bool _has_subdomain_ids
Whether or not subdomain_ids parameter is set.
GeneratedMeshGenerator(const InputParameters ¶meters)
Generates a line, square, or cube mesh with uniformly spaced or biased elements.
std::unique_ptr< T_DEST, T_DELETER > dynamic_pointer_cast(std::unique_ptr< T_SRC, T_DELETER > &src)
These are reworked from https://stackoverflow.com/a/11003103.
const boundary_id_type _boundary_id_offset
offset that is added to the boundary IDs
virtual const std::string & name() const
Get the name of the class.
auto max(const L &left, const R &right)
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
static InputParameters validParams()
registerMooseObject("MooseApp", GeneratedMeshGenerator)
const T & getParam(const std::string &name) const
Retrieve a parameter for the object.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
unsigned int & _nx
Number of elements in x, y, z direction.
Real & _xmin
The min/max values for x,y,z component.
void paramError(const std::string ¶m, Args... args) const
Emits an error prefixed with the file and line number of the given param (from the input file) along ...
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const bool _gauss_lobatto_grid
All of the libmesh build_line/square/cube routines support an option to grade the mesh into the bound...
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
std::unique_ptr< MeshBase > generate() override
Generate / modify the mesh.
MooseEnum _dim
The dimension of the mesh.
const std::string _boundary_name_prefix
prefix string for the boundary names
std::unique_ptr< MeshBase > buildMeshBaseObject(unsigned int dim=libMesh::invalid_uint)
Build a MeshBase object whose underlying type will be determined by the Mesh input file block...
MeshGenerators are objects that can modify or add to an existing mesh.
void ErrorVector unsigned int