Line data Source code
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 : #pragma once
11 :
12 : #include "MeshGenerator.h"
13 : #include "MooseEnum.h"
14 : #include "CastUniquePointer.h"
15 : #include "MooseMeshUtils.h"
16 : #include "libmesh/replicated_mesh.h"
17 : #include "libmesh/mesh_modification.h"
18 : #include "libmesh/face_quad4.h"
19 : #include "libmesh/face_quad8.h"
20 : #include "libmesh/face_quad9.h"
21 : #include "libmesh/face_tri3.h"
22 : #include "libmesh/face_tri6.h"
23 : #include "libmesh/face_tri7.h"
24 : #include "libmesh/serial_mesh.h"
25 : #include "libmesh/boundary_info.h"
26 : #include "libmesh/utility.h"
27 :
28 : /**
29 : * A base class that contains common members for Reactor module mesh generators.
30 : */
31 0 : class PolygonMeshGeneratorBase : public MeshGenerator
32 : {
33 : public:
34 : static InputParameters validParams();
35 :
36 : PolygonMeshGeneratorBase(const InputParameters & parameters);
37 :
38 : virtual std::unique_ptr<MeshBase> generate() override;
39 :
40 : /// An enum class for style of input polygon size
41 : enum class PolygonSizeStyle
42 : {
43 : apothem,
44 : radius
45 : };
46 :
47 : enum MESH_TYPE
48 : {
49 : CORNER_MESH = 1,
50 : BOUNDARY_MESH = 2,
51 : INNER_MESH = 3
52 : };
53 :
54 : enum RETURN_TYPE
55 : {
56 : ANGLE_DEGREE = 1,
57 : ANGLE_TANGENT = 2
58 : };
59 :
60 : enum INTRISIC_SUBDOMAIN_ID : subdomain_id_type
61 : {
62 : PERIPHERAL_ID_SHIFT = 1000,
63 : TRANSITION_LAYER_DEFAULT = 10000
64 : };
65 : // boundary_id_type is short int so max is 32,767
66 : enum INTRINSIC_SIDESET_ID : boundary_id_type
67 : {
68 : OUTER_SIDESET_ID = 10000,
69 : OUTER_SIDESET_ID_ALT = 15000,
70 : SLICE_BEGIN = 30000,
71 : SLICE_END = 31000,
72 : SLICE_ALT = 30500
73 : };
74 :
75 : enum INTRINSIC_NUM_SIDES
76 : {
77 : HEXAGON_NUM_SIDES = 6,
78 : SQUARE_NUM_SIDES = 4
79 : };
80 :
81 : enum class TRI_ELEM_TYPE
82 : {
83 : TRI3,
84 : TRI6,
85 : TRI7
86 : };
87 :
88 : enum class QUAD_ELEM_TYPE
89 : {
90 : QUAD4,
91 : QUAD8,
92 : QUAD9
93 : };
94 :
95 : /// Contains multiple blocks's boundary layer related parameters
96 : struct multiBdryLayerParams
97 : {
98 : std::vector<Real> widths;
99 : std::vector<Real> fractions;
100 : std::vector<unsigned int> intervals;
101 : std::vector<Real> biases;
102 : };
103 : /// Contains a single block's boundary layer related parameters
104 : struct singleBdryLayerParams
105 : {
106 : Real width;
107 : Real fraction;
108 : unsigned int intervals;
109 : Real bias;
110 : };
111 :
112 : protected:
113 : /**
114 : * Creates a mesh of a slice that corresponds to a single side of the polygon to be generated.
115 : * @param ring_radii radii of the ring regions
116 : * @param ring_layers numbers of radial intervals of the ring regions
117 : * @param ring_radial_biases values used for radial meshing biasing in ring regions
118 : * @param ring_inner_boundary_layer_params widths, radial fractions, radial sectors, and growth
119 : * factors of the inner boundary layer of the ring regions
120 : * @param ring_outer_boundary_layer_params widths, radial fractions, radial sectors, and growth
121 : * factors of the outer boundary layer of the ring regions
122 : * @param ducts_center_dist distance parameters of the duct regions
123 : * @param ducts_layers numbers of radial intervals of the duct regions
124 : * @param duct_radial_biases values used for radial meshing biasing in duct regions
125 : * @param duct_inner_boundary_layer_params widths, radial fractions, radial sectors, and growth
126 : * factors of the inner boundary layer of the duct regions
127 : * @param duct_outer_boundary_layer_params widths, radial fractions, radial sectors, and growth
128 : * factors of the outer boundary layer of the duct regions
129 : * @param pitch twice the distance from the ring center vertex to the side defined by the other
130 : * vertices
131 : * @param num_sectors_per_side number of azimuthal intervals
132 : * @param background_intervals number of radial intervals of the background region
133 : * @param background_radial_bias value used for radial meshing biasing in background region
134 : * @param background_inner_boundary_layer_params width, radial sectors, and growth factor of the
135 : * inner boundary layer of the background region
136 : * @param background_outer_boundary_layer_params width, radial sectors, and growth factor of the
137 : * outer boundary layer of the background region
138 : * @param node_id_background_meta pointer to the first node's id of the background region
139 : * @param side_number number of sides of the polygon
140 : * @param side_index index of the polygon side
141 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
142 : * adaptive boundary matching
143 : * @param block_id_shift shift of the subdomain ids generated by this function
144 : * @param quad_center_elements whether the central region contrains quad elements or not
145 : * @param center_quad_factor A fractional radius factor used to determine the radial positions of
146 : * transition nodes in the center region meshed by quad elements (default is 1.0 - 1.0/div_num)
147 : * @param create_inward_interface_boundaries whether inward interface boundary sidesets are
148 : * created
149 : * @param create_outward_interface_boundaries whether outward interface boundary sidesets are
150 : * created
151 : * @param boundary_id_shift shift of the interface boundary ids
152 : * @param generate_side_specific_boundaries whether the side-specific external boundaries are
153 : * generated or not
154 : * @param tri_elem_type type of the triangular elements to be generated
155 : * @param quad_elem_type type of the quadrilateral elements to be generated
156 : * @return a mesh of a polygon slice
157 : */
158 : std::unique_ptr<ReplicatedMesh>
159 : buildSimpleSlice(std::vector<Real> ring_radii,
160 : const std::vector<unsigned int> ring_layers,
161 : const std::vector<Real> ring_radial_biases,
162 : const multiBdryLayerParams & ring_inner_boundary_layer_params,
163 : const multiBdryLayerParams & ring_outer_boundary_layer_params,
164 : std::vector<Real> ducts_center_dist,
165 : const std::vector<unsigned int> ducts_layers,
166 : const std::vector<Real> duct_radial_biases,
167 : const multiBdryLayerParams & duct_inner_boundary_layer_params,
168 : const multiBdryLayerParams & duct_outer_boundary_layer_params,
169 : const Real pitch,
170 : const unsigned int num_sectors_per_side,
171 : const unsigned int background_intervals,
172 : const Real background_radial_bias,
173 : const singleBdryLayerParams & background_inner_boundary_layer_params,
174 : const singleBdryLayerParams & background_outer_boundary_layer_params,
175 : dof_id_type & node_id_background_meta,
176 : const unsigned int side_number,
177 : const unsigned int side_index,
178 : const std::vector<Real> azimuthal_tangent = std::vector<Real>(),
179 : const subdomain_id_type block_id_shift = 0,
180 : const bool quad_center_elements = false,
181 : const Real center_quad_factor = 0.0,
182 : const bool create_inward_interface_boundaries = false,
183 : const bool create_outward_interface_boundaries = true,
184 : const boundary_id_type boundary_id_shift = 0,
185 : const bool generate_side_specific_boundaries = true,
186 : const TRI_ELEM_TYPE tri_elem_type = TRI_ELEM_TYPE::TRI3,
187 : const QUAD_ELEM_TYPE quad_elem_type = QUAD_ELEM_TYPE::QUAD4);
188 :
189 : /**
190 : * Creates a mesh of a general polygon slice with a triangular shape and circular regions on one
191 : * of its vertex.
192 : * @param ring_radii radii of the ring regions
193 : * @param ring_layers numbers of radial intervals of the ring regions
194 : * @param ring_radial_biases values used for radial meshing biasing in ring regions
195 : * @param ring_inner_boundary_layer_params widths, radial fractions, radial sectors, and growth
196 : * factors of the inner boundary layer of the ring regions
197 : * @param ring_outer_boundary_layer_params widths, radial fractions, radial sectors, and growth
198 : * factors of the outer boundary layer of the ring regions
199 : * @param ducts_center_dist distance parameters of the duct regions
200 : * @param ducts_layers numbers of radial intervals of the duct regions
201 : * @param duct_radial_biases values used for radial meshing biasing in duct regions
202 : * @param duct_inner_boundary_layer_params widths, radial fractions, radial sectors, and growth
203 : * factors of the inner boundary layer of the duct regions
204 : * @param duct_outer_boundary_layer_params widths, radial fractions, radial sectors, and growth
205 : * factors of the outer boundary layer of the duct regions
206 : * @param primary_side_length length of the first side (i.e., the side that is parallel to y-axis
207 : * when rotation_angle is zero) that involves the ring center vertex
208 : * @param secondary_side_length length of the second side (obtained by clockwise rotating the fist
209 : * side by azimuthal_angle) that involves the ring center vertex
210 : * @param num_sectors_per_side number of azimuthal intervals
211 : * @param background_intervals number of radial intervals of the background region
212 : * @param background_radial_bias value used for radial meshing biasing in background region
213 : * @param background_inner_boundary_layer_params width, radial sectors, and growth factor of the
214 : * inner boundary layer of the background region
215 : * @param background_outer_boundary_layer_params width, radial sectors, and growth factor of the
216 : * outer boundary layer of the background region
217 : * @param node_id_background_meta pointer to the first node's id of the background region
218 : * @param azimuthal_angle the angle defined by the primary and secondary sides
219 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
220 : * adaptive boundary matching
221 : * @param side_index index of the polygon side
222 : * @param quad_center_elements whether the central region contains quad elements or not
223 : * @param center_quad_factor A fractional radius factor used to determine the radial positions of
224 : * transition nodes in the center region meshed by quad elements (default is 1.0 - 1.0/div_num)
225 : * @param rotation_angle azimuthal angle of the primary side
226 : * @param generate_side_specific_boundaries whether the side-specific external boundaries are
227 : * generated or not
228 : * @return a mesh of a general slice
229 : */
230 : std::unique_ptr<ReplicatedMesh>
231 : buildGeneralSlice(std::vector<Real> ring_radii,
232 : const std::vector<unsigned int> ring_layers,
233 : const std::vector<Real> ring_radial_biases,
234 : const multiBdryLayerParams & ring_inner_boundary_layer_params,
235 : const multiBdryLayerParams & ring_outer_boundary_layer_params,
236 : std::vector<Real> ducts_center_dist,
237 : const std::vector<unsigned int> ducts_layers,
238 : const std::vector<Real> duct_radial_biases,
239 : const multiBdryLayerParams & duct_inner_boundary_layer_params,
240 : const multiBdryLayerParams & duct_outer_boundary_layer_params,
241 : const Real primary_side_length,
242 : const Real secondary_side_length,
243 : const unsigned int num_sectors_per_side,
244 : const unsigned int background_intervals,
245 : const Real background_radial_bias,
246 : const singleBdryLayerParams & background_inner_boundary_layer_params,
247 : const singleBdryLayerParams & background_outer_boundary_layer_params,
248 : dof_id_type & node_id_background_meta,
249 : const Real azimuthal_angle,
250 : const std::vector<Real> azimuthal_tangent,
251 : const unsigned int side_index,
252 : const bool quad_center_elements,
253 : const Real center_quad_factor,
254 : const Real rotation_angle,
255 : const bool generate_side_specific_boundaries = true);
256 :
257 : /**
258 : * Generates a mesh of a polygon slice, which is the foundation of both buildGeneralSlice and
259 : * buildSimpleSlice.
260 : * @param ring_radii radii of the ring regions
261 : * @param ring_layers numbers of radial intervals of the ring regions
262 : * @param ring_radial_biases values used for radial meshing biasing in ring regions
263 : * @param ring_inner_boundary_layer_params widths, radial fractions, radial sectors, and growth
264 : * factors of the inner boundary layer of the ring regions
265 : * @param ring_outer_boundary_layer_params widths, radial fractions, radial sectors, and growth
266 : * factors of the outer boundary layer of the ring regions
267 : * @param ducts_center_dist distance parameters of the duct regions
268 : * @param ducts_layers numbers of radial intervals of the duct regions
269 : * @param duct_radial_biases values used for radial meshing biasing in duct regions
270 : * @param duct_inner_boundary_layer_params widths, radial fractions, radial sectors, and growth
271 : * factors of the inner boundary layer of the duct regions
272 : * @param duct_outer_boundary_layer_params widths, radial fractions, radial sectors, and growth
273 : * factors of the outer boundary layer of the duct regions
274 : * @param pitch twice of the length of the first side times cosine of the azimuthal angle
275 : * @param num_sectors_per_side number of azimuthal intervals
276 : * @param background_intervals number of radial intervals of the background region
277 : * @param background_radial_bias value used for radial meshing biasing in background region
278 : * @param background_inner_boundary_layer_params width, radial sectors, and growth factor of the
279 : * inner boundary layer of the background region
280 : * @param background_outer_boundary_layer_params width, radial sectors, and growth factor of the
281 : * outer boundary layer of the background region
282 : * @param node_id_background_meta pointer to the first node's id of the background region
283 : * @param virtual_side_number 360.0 over the azimuthal angle of the slice (happens to be number of
284 : * sides of the polygon if a regular polygon is to be generated)
285 : * @param side_index index of the polygon side
286 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
287 : * adaptive boundary matching
288 : * @param block_id_shift shift of the subdomain ids generated by this function
289 : * @param quad_center_elements whether the central region contrains quad elements or not
290 : * @param center_quad_factor A fractional radius factor used to determine the radial positions of
291 : * transition nodes in the center region meshed by quad elements (default is 1.0 - 1.0/div_num)
292 : * @param create_inward_interface_boundaries whether inward interface boundary sidesets are
293 : * created
294 : * @param create_outward_interface_boundaries whether outward interface boundary sidesets are
295 : * created
296 : * @param boundary_id_shift shift of the interface boundary ids
297 : * @param pitch_scale_factor the ratio between the secondary side length to the primary side
298 : * length.
299 : * @param generate_side_specific_boundaries whether the side-specific external boundaries are
300 : * generated or not
301 : * @param tri_elem_type type of the triangular elements to be generated
302 : * @param quad_elem_type type of the quadrilateral elements to be generated
303 : * @return a mesh of a slice
304 : */
305 : std::unique_ptr<ReplicatedMesh>
306 : buildSlice(std::vector<Real> ring_radii,
307 : const std::vector<unsigned int> ring_layers,
308 : const std::vector<Real> ring_radial_biases,
309 : const multiBdryLayerParams & ring_inner_boundary_layer_params,
310 : const multiBdryLayerParams & ring_outer_boundary_layer_params,
311 : std::vector<Real> ducts_center_dist,
312 : const std::vector<unsigned int> ducts_layers,
313 : const std::vector<Real> duct_radial_biases,
314 : const multiBdryLayerParams & duct_inner_boundary_layer_params,
315 : const multiBdryLayerParams & duct_outer_boundary_layer_params,
316 : const Real pitch,
317 : const unsigned int num_sectors_per_side,
318 : const unsigned int background_intervals,
319 : const Real background_radial_bias,
320 : const singleBdryLayerParams & background_inner_boundary_layer_params,
321 : const singleBdryLayerParams & background_outer_boundary_layer_params,
322 : dof_id_type & node_id_background_meta,
323 : const Real virtual_side_number,
324 : const unsigned int side_index,
325 : const std::vector<Real> azimuthal_tangent = std::vector<Real>(),
326 : const subdomain_id_type block_id_shift = 0,
327 : const bool quad_center_elements = false,
328 : const Real center_quad_factor = 0.0,
329 : const bool create_inward_interface_boundaries = false,
330 : const bool create_outward_interface_boundaries = true,
331 : const boundary_id_type boundary_id_shift = 0,
332 : const Real pitch_scale_factor = 1.0,
333 : const bool generate_side_specific_boundaries = true,
334 : const TRI_ELEM_TYPE tri_elem_type = TRI_ELEM_TYPE::TRI3,
335 : const QUAD_ELEM_TYPE quad_elem_type = QUAD_ELEM_TYPE::QUAD4);
336 :
337 : /**
338 : * Creates nodes of the very central mesh layer of the polygon for quad central elements.
339 : * @param mesh input mesh to add the nodes onto
340 : * @param virtual_side_number virtual number of sides of the polygon (360/slice_azimuthal)
341 : * @param div_num division number of the central mesh layer
342 : * @param ring_radii_0 radius of the central mesh layer
343 : * @param nodes pointer to the mesh's nodes
344 : * @param nodes vector that contains the nodes with basic geometry information
345 : */
346 : void centerNodes(ReplicatedMesh & mesh,
347 : const Real virtual_side_number,
348 : const unsigned int div_num,
349 : const Real ring_radii_0,
350 : std::vector<std::vector<Node *>> & nodes) const;
351 :
352 : /**
353 : * Creates nodes for the ring-geometry region of a single slice.
354 : * @param mesh input mesh to add the nodes onto
355 : * @param ring_radii radii of the ring regions
356 : * @param ring_layers numbers of radial intervals of the ring regions
357 : * @param biased_terms normalized spacing values used for radial meshing biasing in ring regions
358 : * @param num_sectors_per_side number of azimuthal intervals
359 : * @param corner_p[2][2] array contains the coordinates of the corner positions
360 : * @param corner_to_corner diameter of the circumscribed circle of the polygon
361 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
362 : * adaptive boundary matching
363 : */
364 : void ringNodes(ReplicatedMesh & mesh,
365 : const std::vector<Real> ring_radii,
366 : const std::vector<unsigned int> ring_layers,
367 : const std::vector<std::vector<Real>> biased_terms,
368 : const unsigned int num_sectors_per_side,
369 : const Real corner_p[2][2],
370 : const Real corner_to_corner,
371 : const std::vector<Real> azimuthal_tangent = std::vector<Real>()) const;
372 :
373 : /**
374 : * Creates nodes for the ring-to-polygon transition region (i.e., background) of a single slice.
375 : * @param mesh input mesh to add the nodes onto
376 : * @param num_sectors_per_side number of azimuthal intervals
377 : * @param background_intervals number of radial intervals of the background region
378 : * @param biased_terms normalized spacing values used for radial meshing biasing in background
379 : * region
380 : * @param background_corner_distance center to duct (innermost duct) corner distance
381 : * @param background_corner_radial_interval_length radial interval distance
382 : * @param corner_p[2][2] array contains the coordinates of the corner positions
383 : * @param corner_to_corner diameter of the circumscribed circle of the polygon
384 : * @param background_in radius of the inner boundary of the background region
385 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
386 : * adaptive boundary matching
387 : * @return a mesh with background region nodes created
388 : */
389 : void backgroundNodes(ReplicatedMesh & mesh,
390 : const unsigned int num_sectors_per_side,
391 : const unsigned int background_intervals,
392 : const std::vector<Real> biased_terms,
393 : const Real background_corner_distance,
394 : const Real background_corner_radial_interval_length,
395 : const Real corner_p[2][2],
396 : const Real corner_to_corner,
397 : const Real background_in,
398 : const std::vector<Real> azimuthal_tangent = std::vector<Real>()) const;
399 :
400 : /**
401 : * Creates nodes for the duct-geometry region of a single slice.
402 : * @param mesh input mesh to add the nodes onto
403 : * @param ducts_center_dist distance parameters of the duct regions
404 : * @param ducts_layers numbers of radial intervals of the duct regions
405 : * @param biased_terms normalized spacing values used for radial meshing biasing in duct region
406 : * @param num_sectors_per_side number of azimuthal intervals
407 : * @param corner_p[2][2] array contains the coordinates of the corner positions
408 : * @param corner_to_corner diameter of the circumscribed circle of the polygon
409 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
410 : * adaptive boundary matching
411 : */
412 : void ductNodes(ReplicatedMesh & mesh,
413 : std::vector<Real> * const ducts_center_dist,
414 : const std::vector<unsigned int> ducts_layers,
415 : const std::vector<std::vector<Real>> biased_terms,
416 : const unsigned int num_sectors_per_side,
417 : const Real corner_p[2][2],
418 : const Real corner_to_corner,
419 : const std::vector<Real> azimuthal_tangent = std::vector<Real>()) const;
420 :
421 : /**
422 : * Defines quad elements in the very central region of the polygon.
423 : * @param mesh input mesh to create the elements onto
424 : * @param div_num division number of the central mesh layer
425 : * @param block_id_shift shift of the subdomain ids generated by this function
426 : * @param create_outward_interface_boundaries whether outward interface boundary sidesets are
427 : * created
428 : * @param boundary_id_shift shift of the interface boundary ids
429 : * @param id_array pointer to a vector that contains the node_ids with basic geometry information
430 : * @param assign_external_boundary whether the external boundary ids are assigned
431 : * @param side_index index of the polygon side (only used if external boundary ids are assigned)
432 : * @param generate_side_specific_boundaries whether the side-specific external boundaries are
433 : * generated or not
434 : * @param quad_elem_type type of the quadrilateral elements to be generated
435 : */
436 : void cenQuadElemDef(ReplicatedMesh & mesh,
437 : const unsigned int div_num,
438 : const subdomain_id_type block_id_shift,
439 : const bool create_outward_interface_boundaries,
440 : const boundary_id_type boundary_id_shift,
441 : std::vector<std::vector<Node *>> & nodes,
442 : const bool assign_external_boundary = false,
443 : const unsigned int side_index = 0,
444 : const bool generate_side_specific_boundaries = true,
445 : const QUAD_ELEM_TYPE quad_elem_type = QUAD_ELEM_TYPE::QUAD4) const;
446 :
447 : /**
448 : * Defines triangular elements in the very central region of the polygon.
449 : * @param mesh input mesh to create the elements onto
450 : * @param num_sectors_per_side number of azimuthal intervals
451 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
452 : * adaptive boundary matching
453 : * @param block_id_shift shift of the subdomain ids generated by this function
454 : * @param create_outward_interface_boundaries whether outward interface boundary sidesets are
455 : * created
456 : * @param boundary_id_shift shift of the interface boundary ids
457 : * @param assign_external_boundary whether the external boundary ids are assigned
458 : * @param side_index index of the polygon side (only used if external boundary ids are assigned)
459 : * @param generate_side_specific_boundaries whether the side-specific external boundaries are
460 : * generated or not
461 : * @param tri_elem_type type of the triangular elements to be generated
462 : */
463 : void cenTriElemDef(ReplicatedMesh & mesh,
464 : const unsigned int num_sectors_per_side,
465 : const std::vector<Real> azimuthal_tangent = std::vector<Real>(),
466 : const subdomain_id_type block_id_shift = 0,
467 : const bool create_outward_interface_boundaries = true,
468 : const boundary_id_type boundary_id_shift = 0,
469 : const bool assign_external_boundary = false,
470 : const unsigned int side_index = 0,
471 : const bool generate_side_specific_boundaries = true,
472 : const TRI_ELEM_TYPE tri_elem_type = TRI_ELEM_TYPE::TRI3) const;
473 :
474 : /**
475 : * Defines general quad elements for the polygon.
476 : * @param mesh input mesh to create the elements onto
477 : * @param num_sectors_per_side number of azimuthal intervals
478 : * @param subdomain_rings numbers of radial intervals of all involved subdomain layers
479 : * @param side_index index of the polygon side
480 : * @param azimuthal_tangent vector of tangent values of the azimuthal angles as reference for
481 : * adaptive boundary matching
482 : * @param block_id_shift shift of the subdomain ids generated by this function
483 : * @param nodeid_shift shift of the node_ids of these elements
484 : * @param create_inward_interface_boundaries whether inward interface boundary sidesets are
485 : * created
486 : * @param create_outward_interface_boundaries whether outward interface boundary sidesets are
487 : * created
488 : * @param boundary_id_shift shift of the interface boundary ids
489 : * @param generate_side_specific_boundaries whether the side-specific external boundaries are
490 : * generated or not
491 : * @param quad_elem_type type of the quadrilateral elements to be generated
492 : */
493 : void quadElemDef(ReplicatedMesh & mesh,
494 : const unsigned int num_sectors_per_side,
495 : const std::vector<unsigned int> subdomain_rings,
496 : const unsigned int side_index,
497 : const std::vector<Real> azimuthal_tangent = std::vector<Real>(),
498 : const subdomain_id_type block_id_shift = 0,
499 : const dof_id_type nodeid_shift = 0,
500 : const bool create_inward_interface_boundaries = false,
501 : const bool create_outward_interface_boundaries = true,
502 : const boundary_id_type boundary_id_shift = 0,
503 : const bool generate_side_specific_boundaries = true,
504 : const QUAD_ELEM_TYPE quad_elem_type = QUAD_ELEM_TYPE::QUAD4) const;
505 :
506 : /**
507 : * Creates peripheral area mesh for the patterned hexagon mesh. Note that the function create the
508 : * peripheral area for each side of the unit hexagon mesh before stitching. An edge unit hexagon
509 : * has two sides that need peripheral areas, whereas a corner unit hexagon has three such sides.
510 : * The positions of the inner and outer boundary nodes are pre-calculated as positions_inner and
511 : * d_positions_outer; This function performs interpolation to generate the mesh grid.
512 : * @param mesh input mesh to create the peripheral area mesh onto
513 : * @param num_sectors_per_side number of azimuthal intervals
514 : * @param peripheral_invervals number of radial intervals of the peripheral region
515 : * @param position_inner key positions of the inner side of the peripheral region
516 : * @param d_position_outer key inremental positions of the outer side of the peripheral region
517 : * @param id_shift shift of subdomain id of the peripheral region
518 : * @param create_inward_interface_boundaries whether inward interface boundary sidesets are
519 : * created
520 : * @param create_outward_interface_boundaries whether outward interface boundary sidesets are
521 : * created
522 : * @param quad_elem_type type of quad element to be created
523 : * @return a mesh with the peripheral region added to a hexagon input mesh
524 : */
525 : std::unique_ptr<ReplicatedMesh>
526 : buildSimplePeripheral(const unsigned int num_sectors_per_side,
527 : const unsigned int peripheral_invervals,
528 : const std::vector<std::pair<Real, Real>> & position_inner,
529 : const std::vector<std::pair<Real, Real>> & d_position_outer,
530 : const subdomain_id_type id_shift,
531 : const QUAD_ELEM_TYPE quad_elem_type,
532 : const bool create_inward_interface_boundaries = false,
533 : const bool create_outward_interface_boundaries = true);
534 :
535 : /**
536 : * Adjusts the mid-edge node locations in boundary regions when using quadratic elements with
537 : * uniform boundary node spacing enabled.
538 : * @param out_mesh mesh to be adjusted.
539 : * @param boundary_quad_elem_type boundary quad element type.
540 : */
541 : void adjustPeripheralQuadraticElements(MeshBase & out_mesh,
542 : const QUAD_ELEM_TYPE boundary_quad_elem_type) const;
543 :
544 : /**
545 : * Calculates the point coordinates of within a parallelogram region using linear interpolation.
546 : * @param pi_1_x x coordinate of the first inner side point (parallelogram vertex)
547 : * @param pi_1_y y coordinate of the first inner side point (parallelogram vertex)
548 : * @param po_1_x x coordinate of the first outer side point (parallelogram vertex)
549 : * @param po_1_y y coordinate of the first outer side point (parallelogram vertex)
550 : * @param pi_2_x x coordinate of the second inner side point (parallelogram vertex)
551 : * @param pi_2_y y coordinate of the second inner side point (parallelogram vertex)
552 : * @param po_2_x x coordinate of the second outer side point (parallelogram vertex)
553 : * @param po_2_y y coordinate of the second outer side point (parallelogram vertex)
554 : * @param i passed loop index 1
555 : * @param j passed loop index 2
556 : * @param num_sectors_per_side number of azimuthal intervals
557 : * @param peripheral_invervals number of radial intervals of the peripheral region
558 : * @return an interpolated position within a parallelogram
559 : */
560 : std::pair<Real, Real> pointInterpolate(const Real pi_1_x,
561 : const Real pi_1_y,
562 : const Real po_1_x,
563 : const Real po_1_y,
564 : const Real pi_2_x,
565 : const Real pi_2_y,
566 : const Real po_2_x,
567 : const Real po_2_y,
568 : const unsigned int i,
569 : const unsigned int j,
570 : const unsigned int num_sectors_per_side,
571 : const unsigned int peripheral_intervals) const;
572 :
573 : /**
574 : * Calculates x and y coordinates after rotating by theta angle.
575 : * @param x x coordinate of the node to be rotated
576 : * @param y y coordinate of the node to be rotated
577 : * @param theta rotation angle
578 : * @return n/a
579 : */
580 : void nodeCoordRotate(Real & x, Real & y, const Real theta) const;
581 :
582 : /**
583 : * Deforms peripheral region when the external side of a polygon assembly of stitched meshes
584 : * cuts off the stitched meshes.
585 : * @param mesh input mesh to be deformed
586 : * @param orientation orientation angle of the input mesh (move the deformation direction to y)
587 : * @param y_max_0 original maximum y position
588 : * @param y_max_n maximum y position after deformation
589 : * @param y_min minimum y position that is affected by the deformation
590 : * @param mesh_type whether the peripheral region is for a corner or a side hexagon mesh.
591 : * @param tols tolerance used to determine the boundary of deformation region
592 : * @param unit_angle unit angle of the geometry, which is 60.0 for hexagonal and 90.0 for square
593 : * @return n/a (input mesh is directly altered)
594 : */
595 : void cutOffPolyDeform(MeshBase & mesh,
596 : const Real orientation,
597 : const Real y_max_0,
598 : const Real y_max_n,
599 : const Real y_min,
600 : const unsigned int mesh_type,
601 : const Real unit_angle = 60.0,
602 : const Real tols = 1E-5) const;
603 :
604 : /**
605 : * Finds the center of a quadrilateral based on four vertices.
606 : * @param p1 vertex 1
607 : * @param p2 vertex 2
608 : * @param p3 vertex 3
609 : * @param p4 vertex 4
610 : * @return the intecept point coordinate x and y
611 : */
612 : std::pair<Real, Real> fourPointIntercept(const std::pair<Real, Real> & p1,
613 : const std::pair<Real, Real> & p2,
614 : const std::pair<Real, Real> & p3,
615 : const std::pair<Real, Real> & p4) const;
616 :
617 : /**
618 : * Collects sorted azimuthal angles of the external boundary.
619 : * @param mesh input mesh whose boundary node azimuthal angles need to be collected
620 : * @param boundary_points reference vector to contain the Points corresponding to the collected
621 : * azimuthal angles
622 : * @param lower_azi lower boundary of the azimuthal angles to be collected
623 : * @param upper_azi upper boundary of the azimuthal angles to be collected
624 : * @param return_type whether angle values or tangent values are returned
625 : * @param num_sides number of sides of the input mesh (only used if return type is ANGLE_TANGENT)
626 : * @param bid id of the boundary of which the nodes' azimuthal angles are collected
627 : * @param calculate_origin whether the mesh origin is calculated based on the centroid position
628 : * @param input_origin_x precalculated mesh origin coordinate x
629 : * @param input_origin_y precalculated mesh origin coordinate y
630 : * @param tol tolerance that the minimum azimuthal angle is
631 : * @return the list of azimuthal angles of all the nodes on the external grain boundary within the
632 : * given range
633 : */
634 : std::vector<Real> azimuthalAnglesCollector(ReplicatedMesh & mesh,
635 : std::vector<Point> & boundary_points,
636 : const Real lower_azi = -30.0,
637 : const Real upper_azi = 30.0,
638 : const unsigned int return_type = ANGLE_TANGENT,
639 : const unsigned int num_sides = 6,
640 : const boundary_id_type bid = OUTER_SIDESET_ID,
641 : const bool calculate_origin = true,
642 : const Real input_origin_x = 0.0,
643 : const Real input_origin_y = 0.0,
644 : const Real tol = 1.0E-10) const;
645 :
646 : /**
647 : * Collects sorted azimuthal angles of the external boundary.
648 : * @param mesh input mesh whose boundary node azimuthal angles need to be collected
649 : * @param lower_azi lower boundary of the azimuthal angles to be collected
650 : * @param upper_azi upper boundary of the azimuthal angles to be collected
651 : * @param return_type whether angle values or tangent values are returned
652 : * @param num_sides number of sides of the input mesh (only used if return type is ANGLE_TANGENT)
653 : * @param bid id of the boundary of which the nodes' azimuthal angles are collected
654 : * @param calculate_origin whether the mesh origin is calculated based on the centroid position
655 : * @param input_origin_x precalculated mesh origin coordinate x
656 : * @param input_origin_y precalculated mesh origin coordinate y
657 : * @param tol tolerence that the minimum azimuthal angle is
658 : * @return the list of azimuthal angles of all the nodes on the external grain boundary within the
659 : * given range
660 : */
661 : std::vector<Real> azimuthalAnglesCollector(ReplicatedMesh & mesh,
662 : const Real lower_azi = -30.0,
663 : const Real upper_azi = 30.0,
664 : const unsigned int return_type = ANGLE_TANGENT,
665 : const unsigned int num_sides = 6,
666 : const boundary_id_type bid = OUTER_SIDESET_ID,
667 : const bool calculate_origin = true,
668 : const Real input_origin_x = 0.0,
669 : const Real input_origin_y = 0.0,
670 : const Real tol = 1.0E-10) const;
671 :
672 : /**
673 : * Creates bias terms for multiple blocks.
674 : * @param radial_biases bias growth factors of the elements within the main regions of the blocks
675 : * @param intervals radial interval numbers of the main regions of the blocks
676 : * @param inner_boundary_layer_params widths, radial fractions, radial sectors, and growth
677 : * factors of the inner boundary layers
678 : * @param outer_boundary_layer_params widths, radial fractions, radial sectors, and growth
679 : * factors of the outer boundary layers
680 : * @return bias list of terms describing the cumulative radial fractions of the nodes within
681 : * multiple blocks
682 : */
683 : std::vector<std::vector<Real>>
684 : biasTermsCalculator(const std::vector<Real> radial_biases,
685 : const std::vector<unsigned int> intervals,
686 : const multiBdryLayerParams inner_boundary_layer_params,
687 : const multiBdryLayerParams outer_boundary_layer_params) const;
688 :
689 : /**
690 : * Creates bias terms for a single block.
691 : * @param radial_bias bias growth factor of the elements within the main region of the block
692 : * @param intervals radial interval number of the main region of the block
693 : * @param inner_boundary_layer_params width, radial fraction, radial sector, and growth
694 : * factor of the inner boundary layer
695 : * @param outer_boundary_layer_params width, radial fraction, radial sector, and growth
696 : * factor of the outer boundary layer
697 : * @return bias terms describing the cumulative radial fractions of the nodes within a single
698 : * block
699 : */
700 : std::vector<Real> biasTermsCalculator(
701 : const Real radial_bias,
702 : const unsigned int intervals,
703 : const singleBdryLayerParams inner_boundary_layer_params = {0.0, 0.0, 0, 1.0},
704 : const singleBdryLayerParams outer_boundary_layer_params = {0.0, 0.0, 0, 1.0}) const;
705 :
706 : /**
707 : * Add InputParameters which are used by ring and sector IDs
708 : * @param params InputParameters to be modified with the added params
709 : */
710 : static void addRingAndSectorIDParams(InputParameters & params);
711 :
712 : /**
713 : * assign sector extra ids to polygon mesh
714 : * @param mesh input mesh where sector extra ids are assigned
715 : * @param id_name sector extra ID name
716 : * @param num_side number of polygon sides
717 : * @param num_sectors_per_side number of sections of each side of the polygon
718 : */
719 : void setSectorExtraIDs(MeshBase & mesh,
720 : const std::string id_name,
721 : const unsigned int num_sides,
722 : const std::vector<unsigned int> num_sectors_per_side);
723 :
724 : /**
725 : * assign ring extra ids to polygon mesh
726 : * @param mesh input mesh where ring extra ids are assigned
727 : * @param id_name ring extra id name
728 : * @param num_sides number of polygon sides
729 : * @param num_sectors_per_side number of sectors of each side of the polygon
730 : * @param ring_intervals number of rings in each circle
731 : * @param ring_wise_id whether ring ids are assigned to each ring or to each block
732 : * @param quad_center_elements whether center elements are quad or triangular
733 : */
734 : void setRingExtraIDs(MeshBase & mesh,
735 : const std::string id_name,
736 : const unsigned int num_sides,
737 : const std::vector<unsigned int> num_sectors_per_side,
738 : const std::vector<unsigned int> ring_intervals,
739 : const bool ring_wise_id,
740 : const bool quad_center_elements);
741 :
742 : /**
743 : * reassign interface boundary IDs on the input mesh by applying the boundary ID shift
744 : * @param mesh input mesh
745 : * @param id_shift ID shift value to be applied
746 : * @param boundary_ids list of boundary IDs to be reassigned
747 : * @param reverse remove boundary ID shift
748 : */
749 : void reassignBoundaryIDs(MeshBase & mesh,
750 : const boundary_id_type id_shift,
751 : const std::set<boundary_id_type> & boundary_ids,
752 : const bool reverse = false);
753 :
754 : /**
755 : * returns a list of interface boundary IDs on the mesh generated by this mesh generator
756 : * @param pattern pattern of cells used in this mesh generator
757 : * @param interface_boundary_id_shift_pattern 2D pattern of shift values applied to the boundary
758 : * IDs inside each pattern cells
759 : * @param boundary_ids list of boundary IDs on the mesh generated by this mesh generator
760 : * @param input_interface_boundary_ids list of interface boundary IDs of the pattern cells
761 : * @param use_interface_boundary_id_shift whether ID shifts are applied to interface boundary IDs
762 : * of the pattern cells
763 : * @param create_interface_boundary_id whether interface boundary IDs are generated by this
764 : * mesh generator
765 : * @param num_extra_layers number of extra layers to define background and duct regions on the
766 : * patterned mesh generated by this mesh generator
767 : */
768 : std::set<boundary_id_type> getInterfaceBoundaryIDs(
769 : const std::vector<std::vector<unsigned int>> & pattern,
770 : const std::vector<std::vector<boundary_id_type>> & interface_boundary_id_shift_pattern,
771 : const std::set<boundary_id_type> & boundary_ids,
772 : const std::vector<std::set<boundary_id_type>> & input_interface_boundary_ids,
773 : const bool use_interface_boundary_id_shift,
774 : const bool create_interface_boundary_id,
775 : const unsigned int num_extra_layers) const;
776 :
777 : /**
778 : * Modifies the input multi boundary layer parameters for node generation, especially for the
779 : * quadratic elements
780 : * @param original_multi_bdry_layer_params original multi boundary layer parameters
781 : * @param order order of the elements
782 : * @return modified multi boundary layer parameters
783 : */
784 : multiBdryLayerParams
785 : modifiedMultiBdryLayerParamsCreator(const multiBdryLayerParams & original_multi_bdry_layer_params,
786 : const unsigned int order) const;
787 :
788 : /**
789 : * Modifies the input single boundary layer parameters for node generation, especially for the
790 : * quadratic elements
791 : * @param original_single_bdry_layer_params original single boundary layer parameters
792 : * @param order order of the elements
793 : * @return modified single boundary layer parameters
794 : */
795 : singleBdryLayerParams modifiedSingleBdryLayerParamsCreator(
796 : const singleBdryLayerParams & original_single_bdry_layer_params,
797 : const unsigned int order) const;
798 :
799 : /**
800 : * Generate a string that contains the detailed metadata information for inconsistent input mesh
801 : * metadata error messages
802 : * @param input_names list of input mesh generator names
803 : * @param metadata_vals list of input mesh metadata values
804 : * @param metadata_name name of the input mesh metadata
805 : * @return a string that contains the detailed metadata information
806 : */
807 : std::string pitchMetaDataErrorGenerator(const std::vector<MeshGeneratorName> & input_names,
808 : const std::vector<Real> & metadata_vals,
809 : const std::string & metadata_name) const;
810 : };
|