Loading [MathJax]/extensions/tex2jax.js
https://mooseframework.inl.gov
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends
GeneratedMesh.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
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");
35  params.addRequiredParam<MooseEnum>("dim", dims, "The dimension of the mesh to be generated");
36 
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");
46  params.addParam<MooseEnum>("elem_type",
47  elem_types,
48  "The type of element from libMesh to "
49  "generate (default: linear element for "
50  "requested dimension)");
51  params.addParam<bool>(
52  "gauss_lobatto_grid",
53  false,
54  "Grade mesh into boundaries according to Gauss-Lobatto quadrature spacing.");
55  params.addRangeCheckedParam<Real>(
56  "bias_x",
57  1.,
58  "bias_x>=0.5 & bias_x<=2",
59  "The amount by which to grow (or shrink) the cells in the x-direction.");
60  params.addRangeCheckedParam<Real>(
61  "bias_y",
62  1.,
63  "bias_y>=0.5 & bias_y<=2",
64  "The amount by which to grow (or shrink) the cells in the y-direction.");
65  params.addRangeCheckedParam<Real>(
66  "bias_z",
67  1.,
68  "bias_z>=0.5 & bias_z<=2",
69  "The amount by which to grow (or shrink) the cells in the z-direction.");
70 
71  params.addParamNamesToGroup("dim", "Required");
72 
73  params.addClassDescription(
74  "Create a line, square, or cube mesh with uniformly spaced or biased elements.");
75  return params;
76 }
77 
79  : MooseMesh(parameters),
80  _dim(getParam<MooseEnum>("dim")),
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)
95 {
96  if (_gauss_lobatto_grid && (_bias_x != 1.0 || _bias_y != 1.0 || _bias_z != 1.0))
97  mooseError("Cannot apply both Gauss-Lobatto mesh grading and biasing at the same time.");
98 
99  // All generated meshes are regular orthogonal meshes - until they get modified ;)
101 }
102 
103 void
105 {
106  MooseMesh::prepared(state);
107 
108  // Fall back on scanning the mesh for coordinates instead of using input parameters for queries
109  if (!state)
110  _dims_may_have_changed = true;
111 }
112 
113 Real
114 GeneratedMesh::getMinInDimension(unsigned int component) const
115 {
117  return MooseMesh::getMinInDimension(component);
118 
119  switch (component)
120  {
121  case 0:
122  return _xmin;
123  case 1:
124  return _dim > 1 ? _ymin : 0;
125  case 2:
126  return _dim > 2 ? _zmin : 0;
127  default:
128  mooseError("Invalid component");
129  }
130 }
131 
132 Real
133 GeneratedMesh::getMaxInDimension(unsigned int component) const
134 {
136  return MooseMesh::getMaxInDimension(component);
137 
138  switch (component)
139  {
140  case 0:
141  return _xmax;
142  case 1:
143  return _dim > 1 ? _ymax : 0;
144  case 2:
145  return _dim > 2 ? _zmax : 0;
146  default:
147  mooseError("Invalid component");
148  }
149 }
150 
151 std::unique_ptr<MooseMesh>
153 {
154  return _app.getFactory().copyConstruct(*this);
155 }
156 
157 void
159 {
160  MooseEnum elem_type_enum = getParam<MooseEnum>("elem_type");
161 
162  if (!isParamValid("elem_type"))
163  {
164  // Switching on MooseEnum
165  switch (_dim)
166  {
167  case 1:
168  elem_type_enum = "EDGE2";
169  break;
170  case 2:
171  elem_type_enum = "QUAD4";
172  break;
173  case 3:
174  elem_type_enum = "HEX8";
175  break;
176  }
177  }
178 
179  ElemType elem_type = Utility::string_to_enum<ElemType>(elem_type_enum);
180 
181  // Switching on MooseEnum
182  switch (_dim)
183  {
184  // The build_XYZ mesh generation functions take an
185  // UnstructuredMesh& as the first argument, hence the dynamic_cast.
186  case 1:
187  MeshTools::Generation::build_line(dynamic_cast<UnstructuredMesh &>(getMesh()),
188  _nx,
189  _xmin,
190  _xmax,
191  elem_type,
193  break;
194  case 2:
195  MeshTools::Generation::build_square(dynamic_cast<UnstructuredMesh &>(getMesh()),
196  _nx,
197  _ny,
198  _xmin,
199  _xmax,
200  _ymin,
201  _ymax,
202  elem_type,
204  break;
205  case 3:
206  MeshTools::Generation::build_cube(dynamic_cast<UnstructuredMesh &>(getMesh()),
207  _nx,
208  _ny,
209  _nz,
210  _xmin,
211  _xmax,
212  _ymin,
213  _ymax,
214  _zmin,
215  _zmax,
216  elem_type,
218  break;
219  }
220 
221  // Apply the bias if any exists
222  if (_bias_x != 1.0 || _bias_y != 1.0 || _bias_z != 1.0)
223  {
224  // Reference to the libmesh mesh
225  MeshBase & mesh = getMesh();
226 
227  // Biases
228  std::array<Real, LIBMESH_DIM> bias = {
229  {_bias_x, _dim > 1 ? _bias_y : 1.0, _dim > 2 ? _bias_z : 1.0}};
230 
231  // "width" of the mesh in each direction
232  std::array<Real, LIBMESH_DIM> width = {
233  {_xmax - _xmin, _dim > 1 ? _ymax - _ymin : 0, _dim > 2 ? _zmax - _zmin : 0}};
234 
235  // Min mesh extent in each direction.
236  std::array<Real, LIBMESH_DIM> mins = {
238 
239  // Number of elements in each direction.
240  std::array<unsigned int, LIBMESH_DIM> nelem = {{_nx, _dim > 1 ? _ny : 1, _dim > 2 ? _nz : 1}};
241 
242  // We will need the biases raised to integer powers in each
243  // direction, so let's pre-compute those...
244  std::array<std::vector<Real>, LIBMESH_DIM> pows;
245  for (unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
246  {
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));
250  }
251 
252  // Loop over the nodes and move them to the desired location
253  for (auto & node_ptr : mesh.node_ptr_range())
254  {
255  Node & node = *node_ptr;
256 
257  for (unsigned int dir = 0; dir < LIBMESH_DIM; ++dir)
258  {
259  if (width[dir] != 0. && bias[dir] != 1.)
260  {
261  // Compute the scaled "index" of the current point. This
262  // will either be close to a whole integer or a whole
263  // integer+0.5 for quadratic nodes.
264  Real float_index = (node(dir) - mins[dir]) * nelem[dir] / width[dir];
265 
266  Real integer_part = 0;
267  Real fractional_part = std::modf(float_index, &integer_part);
268 
269  // Figure out where to move the node...
270  if (std::abs(fractional_part) < TOLERANCE || std::abs(fractional_part - 1.0) < TOLERANCE)
271  {
272  // If the fractional_part ~ 0.0 or 1.0, this is a vertex node, so
273  // we don't need an average.
274  //
275  // Compute the "index" we are at in the current direction. We
276  // round to the nearest integral value to get this instead
277  // of using "integer_part", since that could be off by a
278  // lot (e.g. we want 3.9999 to map to 4.0 instead of 3.0).
279  int index = round(float_index);
280 
281  mooseAssert(index >= static_cast<int>(0) && index < static_cast<int>(pows[dir].size()),
282  "Scaled \"index\" out of range");
283 
284  // Move node to biased location.
285  node(dir) =
286  mins[dir] + width[dir] * (1. - pows[dir][index]) / (1. - pows[dir][nelem[dir]]);
287  }
288  else if (std::abs(fractional_part - 0.5) < TOLERANCE)
289  {
290  // If the fractional_part ~ 0.5, this is a midedge/face
291  // (i.e. quadratic) node. We don't move those with the same
292  // bias as the vertices, instead we put them midway between
293  // their respective vertices.
294  //
295  // Also, since the fractional part is nearly 0.5, we know that
296  // the integer_part will be the index of the vertex to the
297  // left, and integer_part+1 will be the index of the
298  // vertex to the right.
299  node(dir) = mins[dir] +
300  width[dir] *
301  (1. - 0.5 * (pows[dir][integer_part] + pows[dir][integer_part + 1])) /
302  (1. - pows[dir][nelem[dir]]);
303  }
304  else
305  {
306  // We don't yet handle anything higher order than quadratic...
307  mooseError("Unable to bias node at node(", dir, ")=", node(dir));
308  }
309  }
310  }
311  }
312  }
313 }
static InputParameters validParams()
Typical "Moose-style" constructor and copy constructor.
Definition: MooseMesh.C:83
virtual Real getMaxInDimension(unsigned int component) const
Definition: MooseMesh.C:2197
MetaPhysicL::DualNumber< V, D, asd > abs(const MetaPhysicL::DualNumber< V, D, asd > &a)
Definition: EigenADReal.h:42
ElemType
virtual Real getMinInDimension(unsigned int component) const
Returns the min or max of the requested dimension respectively.
Definition: MooseMesh.C:2188
void elem_types(const MeshBase &mesh, std::vector< ElemType > &et)
bool prepared() const
Setter/getter for whether the mesh is prepared.
Definition: MooseMesh.C:3126
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:434
bool isParamValid(const std::string &name) const
Test if the supplied parameter is valid.
Mesh generated from parameters.
Definition: GeneratedMesh.h:17
MeshBase & getMesh()
Accessor for the underlying libMesh Mesh object.
Definition: MooseMesh.C:3437
T round(T x)
Definition: MathUtils.h:77
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:33
Real _xmin
The min/max values for x,y,z component.
Definition: GeneratedMesh.h:43
GeneratedMesh(const InputParameters &parameters)
Definition: GeneratedMesh.C:78
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:817
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 optional 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:1562
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