www.mooseframework.org
GeneratedMesh.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://www.mooseframework.org
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #include "GeneratedMesh.h"
11 
12 #include "MooseApp.h"
13 
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"
20 
21 // C++ includes
22 #include <cmath> // provides round, not std::round (see http://www.cplusplus.com/reference/cmath/round/)
23 #include <array>
24 
26 
29 {
31 
32  MooseEnum elem_types(LIST_GEOM_ELEM); // no default
33 
34  MooseEnum dims("1=1 2 3");
36  "dim", dims, "The dimension of the mesh to be generated"); // Make this parameter required
37 
38  params.addParam<unsigned int>("nx", 1, "Number of elements in the X direction");
39  params.addParam<unsigned int>("ny", 1, "Number of elements in the Y direction");
40  params.addParam<unsigned int>("nz", 1, "Number of elements in the Z direction");
41  params.addParam<Real>("xmin", 0.0, "Lower X Coordinate of the generated mesh");
42  params.addParam<Real>("ymin", 0.0, "Lower Y Coordinate of the generated mesh");
43  params.addParam<Real>("zmin", 0.0, "Lower Z Coordinate of the generated mesh");
44  params.addParam<Real>("xmax", 1.0, "Upper X Coordinate of the generated mesh");
45  params.addParam<Real>("ymax", 1.0, "Upper Y Coordinate of the generated mesh");
46  params.addParam<Real>("zmax", 1.0, "Upper Z Coordinate of the generated mesh");
47  params.addParam<MooseEnum>("elem_type",
48  elem_types,
49  "The type of element from libMesh to "
50  "generate (default: linear element for "
51  "requested dimension)");
52  params.addParam<bool>(
53  "gauss_lobatto_grid",
54  false,
55  "Grade mesh into boundaries according to Gauss-Lobatto quadrature spacing.");
56  params.addRangeCheckedParam<Real>(
57  "bias_x",
58  1.,
59  "bias_x>=0.5 & bias_x<=2",
60  "The amount by which to grow (or shrink) the cells in the x-direction.");
61  params.addRangeCheckedParam<Real>(
62  "bias_y",
63  1.,
64  "bias_y>=0.5 & bias_y<=2",
65  "The amount by which to grow (or shrink) the cells in the y-direction.");
66  params.addRangeCheckedParam<Real>(
67  "bias_z",
68  1.,
69  "bias_z>=0.5 & bias_z<=2",
70  "The amount by which to grow (or shrink) the cells in the z-direction.");
71 
72  params.addParamNamesToGroup("dim", "Required");
73 
74  params.addClassDescription(
75  "Create a line, square, or cube mesh with uniformly spaced or biased elements.");
76  return params;
77 }
78 
80  : MooseMesh(parameters),
81  _dim(getParam<MooseEnum>("dim")),
82  _nx(getParam<unsigned int>("nx")),
83  _ny(getParam<unsigned int>("ny")),
84  _nz(getParam<unsigned int>("nz")),
85  _xmin(getParam<Real>("xmin")),
86  _xmax(getParam<Real>("xmax")),
87  _ymin(getParam<Real>("ymin")),
88  _ymax(getParam<Real>("ymax")),
89  _zmin(getParam<Real>("zmin")),
90  _zmax(getParam<Real>("zmax")),
91  _gauss_lobatto_grid(getParam<bool>("gauss_lobatto_grid")),
92  _bias_x(getParam<Real>("bias_x")),
93  _bias_y(getParam<Real>("bias_y")),
94  _bias_z(getParam<Real>("bias_z")),
95  _dims_may_have_changed(false)
96 {
97  if (_gauss_lobatto_grid && (_bias_x != 1.0 || _bias_y != 1.0 || _bias_z != 1.0))
98  mooseError("Cannot apply both Gauss-Lobatto mesh grading and biasing at the same time.");
99 
100  // All generated meshes are regular orthogonal meshes - until they get modified ;)
102 }
103 
104 void
106 {
107  MooseMesh::prepared(state);
108 
109  // Fall back on scanning the mesh for coordinates instead of using input parameters for queries
110  if (!state)
111  _dims_may_have_changed = true;
112 }
113 
114 Real
115 GeneratedMesh::getMinInDimension(unsigned int component) const
116 {
118  return MooseMesh::getMinInDimension(component);
119 
120  switch (component)
121  {
122  case 0:
123  return _xmin;
124  case 1:
125  return _dim > 1 ? _ymin : 0;
126  case 2:
127  return _dim > 2 ? _zmin : 0;
128  default:
129  mooseError("Invalid component");
130  }
131 }
132 
133 Real
134 GeneratedMesh::getMaxInDimension(unsigned int component) const
135 {
137  return MooseMesh::getMaxInDimension(component);
138 
139  switch (component)
140  {
141  case 0:
142  return _xmax;
143  case 1:
144  return _dim > 1 ? _ymax : 0;
145  case 2:
146  return _dim > 2 ? _zmax : 0;
147  default:
148  mooseError("Invalid component");
149  }
150 }
151 
152 std::unique_ptr<MooseMesh>
154 {
155  return _app.getFactory().copyConstruct(*this);
156 }
157 
158 void
160 {
161  MooseEnum elem_type_enum = getParam<MooseEnum>("elem_type");
162 
163  if (!isParamValid("elem_type"))
164  {
165  // Switching on MooseEnum
166  switch (_dim)
167  {
168  case 1:
169  elem_type_enum = "EDGE2";
170  break;
171  case 2:
172  elem_type_enum = "QUAD4";
173  break;
174  case 3:
175  elem_type_enum = "HEX8";
176  break;
177  }
178  }
179 
180  ElemType elem_type = Utility::string_to_enum<ElemType>(elem_type_enum);
181 
182  // Switching on MooseEnum
183  switch (_dim)
184  {
185  // The build_XYZ mesh generation functions take an
186  // UnstructuredMesh& as the first argument, hence the dynamic_cast.
187  case 1:
188  MeshTools::Generation::build_line(dynamic_cast<UnstructuredMesh &>(getMesh()),
189  _nx,
190  _xmin,
191  _xmax,
192  elem_type,
194  break;
195  case 2:
196  MeshTools::Generation::build_square(dynamic_cast<UnstructuredMesh &>(getMesh()),
197  _nx,
198  _ny,
199  _xmin,
200  _xmax,
201  _ymin,
202  _ymax,
203  elem_type,
205  break;
206  case 3:
207  MeshTools::Generation::build_cube(dynamic_cast<UnstructuredMesh &>(getMesh()),
208  _nx,
209  _ny,
210  _nz,
211  _xmin,
212  _xmax,
213  _ymin,
214  _ymax,
215  _zmin,
216  _zmax,
217  elem_type,
219  break;
220  }
221 
222  // Apply the bias if any exists
223  if (_bias_x != 1.0 || _bias_y != 1.0 || _bias_z != 1.0)
224  {
225  // Reference to the libmesh mesh
226  MeshBase & mesh = getMesh();
227 
228  // Biases
229  std::array<Real, LIBMESH_DIM> bias = {
230  {_bias_x, _dim > 1 ? _bias_y : 1.0, _dim > 2 ? _bias_z : 1.0}};
231 
232  // "width" of the mesh in each direction
233  std::array<Real, LIBMESH_DIM> width = {
234  {_xmax - _xmin, _dim > 1 ? _ymax - _ymin : 0, _dim > 2 ? _zmax - _zmin : 0}};
235 
236  // Min mesh extent in each direction.
237  std::array<Real, LIBMESH_DIM> mins = {
239 
240  // Number of elements in each direction.
241  std::array<unsigned int, LIBMESH_DIM> nelem = {{_nx, _dim > 1 ? _ny : 1, _dim > 2 ? _nz : 1}};
242 
243  // We will need the biases raised to integer powers in each
244  // direction, so let's pre-compute those...
245  std::array<std::vector<Real>, LIBMESH_DIM> pows;
246  for (unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
247  {
248  pows[dir].resize(nelem[dir] + 1);
249  for (unsigned int i = 0; i < pows[dir].size(); ++i)
250  pows[dir][i] = std::pow(bias[dir], static_cast<int>(i));
251  }
252 
253  // Loop over the nodes and move them to the desired location
254  for (auto & node_ptr : mesh.node_ptr_range())
255  {
256  Node & node = *node_ptr;
257 
258  for (unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
259  {
260  if (width[dir] != 0. && bias[dir] != 1.)
261  {
262  // Compute the scaled "index" of the current point. This
263  // will either be close to a whole integer or a whole
264  // integer+0.5 for quadratic nodes.
265  Real float_index = (node(dir) - mins[dir]) * nelem[dir] / width[dir];
266 
267  Real integer_part = 0;
268  Real fractional_part = std::modf(float_index, &integer_part);
269 
270  // Figure out where to move the node...
271  if (std::abs(fractional_part) < TOLERANCE || std::abs(fractional_part - 1.0) < TOLERANCE)
272  {
273  // If the fractional_part ~ 0.0 or 1.0, this is a vertex node, so
274  // we don't need an average.
275  //
276  // Compute the "index" we are at in the current direction. We
277  // round to the nearest integral value to get this instead
278  // of using "integer_part", since that could be off by a
279  // lot (e.g. we want 3.9999 to map to 4.0 instead of 3.0).
280  int index = round(float_index);
281 
282  mooseAssert(index >= static_cast<int>(0) && index < static_cast<int>(pows[dir].size()),
283  "Scaled \"index\" out of range");
284 
285  // Move node to biased location.
286  node(dir) =
287  mins[dir] + width[dir] * (1. - pows[dir][index]) / (1. - pows[dir][nelem[dir]]);
288  }
289  else if (std::abs(fractional_part - 0.5) < TOLERANCE)
290  {
291  // If the fractional_part ~ 0.5, this is a midedge/face
292  // (i.e. quadratic) node. We don't move those with the same
293  // bias as the vertices, instead we put them midway between
294  // their respective vertices.
295  //
296  // Also, since the fractional part is nearly 0.5, we know that
297  // the integer_part will be the index of the vertex to the
298  // left, and integer_part+1 will be the index of the
299  // vertex to the right.
300  node(dir) = mins[dir] +
301  width[dir] *
302  (1. - 0.5 * (pows[dir][integer_part] + pows[dir][integer_part + 1])) /
303  (1. - pows[dir][nelem[dir]]);
304  }
305  else
306  {
307  // We don't yet handle anything higher order than quadratic...
308  mooseError("Unable to bias node at node(", dir, ")=", node(dir));
309  }
310  }
311  }
312  }
313  }
314 }
static InputParameters validParams()
Typical "Moose-style" constructor and copy constructor.
Definition: MooseMesh.C:79
virtual Real getMaxInDimension(unsigned int component) const
Definition: MooseMesh.C:1965
ElemType
virtual Real getMinInDimension(unsigned int component) const
Returns the min or max of the requested dimension respectively.
Definition: MooseMesh.C:1956
void elem_types(const MeshBase &mesh, std::vector< ElemType > &et)
bool prepared() const
Setter/getter for whether the mesh is prepared.
Definition: MooseMesh.C:2888
const std::string LIST_GEOM_ELEM
Definition: MooseMesh.h:58
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...
MeshBase & mesh
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
MooseEnum _dim
The dimension of the mesh.
Definition: GeneratedMesh.h:37
virtual Real getMaxInDimension(unsigned int component) const override
bool _dims_may_have_changed
Boolean to indicate that dimensions may have changed.
Definition: GeneratedMesh.h:60
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.
Definition: Factory.h:310
void addRequiredParam(const std::string &name, const std::string &doc_string)
This method adds a parameter and documentation string to the InputParameters object that will be extr...
virtual void buildMesh() override
Must be overridden by child classes.
Factory & getFactory()
Retrieve a writable reference to the Factory associated with this App.
Definition: MooseApp.h:396
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
ADRealEigenVector< T, D, asd > abs(const ADRealEigenVector< T, D, asd > &)
Mesh generated from parameters.
Definition: GeneratedMesh.h:17
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:3199
T round(T x)
Definition: MathUtils.h:76
MooseMesh wraps a libMesh::Mesh object and enhances its capabilities by caching additional data and s...
Definition: MooseMesh.h:88
unsigned int _ny
Definition: GeneratedMesh.h:40
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:31
Real _xmin
The min/max values for x,y,z component.
Definition: GeneratedMesh.h:43
GeneratedMesh(const InputParameters &parameters)
Definition: GeneratedMesh.C:79
MooseApp & _app
The MOOSE application this is associated with.
Definition: MooseBase.h:84
registerMooseObject("MooseApp", GeneratedMesh)
unsigned int _nz
Definition: GeneratedMesh.h:40
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...
Definition: GeneratedMesh.h:50
virtual const Node & node(const dof_id_type i) const
Various accessors (pointers/references) for Node "i".
Definition: MooseMesh.C:623
void mooseError(Args &&... args) const
Emits an error prefixed with object name and type.
void addClassDescription(const std::string &doc_string)
This method adds a description of the class that will be displayed in the input file syntax dump...
unsigned int _nx
Number of elements in x, y, z direction.
Definition: GeneratedMesh.h:40
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object...
void addRangeCheckedParam(const std::string &name, const T &value, const std::string &parsed_function, const std::string &doc_string)
Real _bias_x
The amount by which to bias the cells in the x,y,z directions.
Definition: GeneratedMesh.h:57
MooseUnits pow(const MooseUnits &, int)
Definition: Units.C:537
void ErrorVector unsigned int
bool _regular_orthogonal_mesh
Boolean indicating whether this mesh was detected to be regular and orthogonal.
Definition: MooseMesh.h:1533
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...
static InputParameters validParams()
Definition: GeneratedMesh.C:28