https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Static Public Attributes | Protected Member Functions | Static Protected Member Functions | Protected Attributes | List of all members
SCMDetailedTriSubChannelMeshGenerator Class Reference

Mesh generator that builds a 3D mesh representing triangular subchannels. More...

#include <SCMDetailedTriSubChannelMeshGenerator.h>

Inheritance diagram for SCMDetailedTriSubChannelMeshGenerator:
[legend]

Public Types

typedef DataFileName DataFileParameterType
 

Public Member Functions

 SCMDetailedTriSubChannelMeshGenerator (const InputParameters &parameters)
 
virtual std::unique_ptr< MeshBase > generate () override
 
std::unique_ptr< MeshBase > generateInternal ()
 
const std::set< MeshGeneratorName > & getRequestedMeshGenerators () const
 
const std::set< MeshGeneratorName > & getRequestedMeshGeneratorsForSub () const
 
void addParentMeshGenerator (const MeshGenerator &mg, const AddParentChildKey)
 
void addChildMeshGenerator (const MeshGenerator &mg, const AddParentChildKey)
 
const std::set< const MeshGenerator *, Comparator > & getParentMeshGenerators () const
 
const std::set< const MeshGenerator *, Comparator > & getChildMeshGenerators () const
 
const std::set< const MeshGenerator *, Comparator > & getSubMeshGenerators () const
 
bool isParentMeshGenerator (const MeshGeneratorName &name, const bool direct=true) const
 
bool isChildMeshGenerator (const MeshGeneratorName &name, const bool direct=true) const
 
bool isNullMeshName (const MeshGeneratorName &name) const
 
bool hasSaveMesh () const
 
bool hasOutput () const
 
const std::string & getSavedMeshName () const
 
bool hasGenerateData () const
 
bool isDataOnly () const
 
virtual bool enabled () const
 
std::shared_ptr< MooseObjectgetSharedPtr ()
 
std::shared_ptr< const MooseObjectgetSharedPtr () const
 
MooseAppgetMooseApp () const
 
const std::string & type () const
 
virtual const std::string & name () const
 
std::string typeAndName () const
 
std::string errorPrefix (const std::string &error_type) const
 
void callMooseError (std::string msg, const bool with_prefix) const
 
MooseObjectParameterName uniqueParameterName (const std::string &parameter_name) const
 
const InputParametersparameters () const
 
MooseObjectName uniqueName () const
 
const T & getParam (const std::string &name) const
 
std::vector< std::pair< T1, T2 > > getParam (const std::string &param1, const std::string &param2) const
 
const T * queryParam (const std::string &name) const
 
const T & getRenamedParam (const std::string &old_name, const std::string &new_name) const
 
getCheckedPointerParam (const std::string &name, const std::string &error_string="") const
 
bool isParamValid (const std::string &name) const
 
bool isParamSetByUser (const std::string &nm) const
 
void paramError (const std::string &param, Args... args) const
 
void paramWarning (const std::string &param, Args... args) const
 
void paramInfo (const std::string &param, Args... args) const
 
void connectControllableParams (const std::string &parameter, const std::string &object_type, const std::string &object_name, const std::string &object_parameter) const
 
void mooseError (Args &&... args) const
 
void mooseErrorNonPrefixed (Args &&... args) const
 
void mooseDocumentedError (const std::string &repo_name, const unsigned int issue_num, Args &&... args) const
 
void mooseWarning (Args &&... args) const
 
void mooseWarningNonPrefixed (Args &&... args) const
 
void mooseDeprecated (Args &&... args) const
 
void mooseInfo (Args &&... args) const
 
std::string getDataFileName (const std::string &param) const
 
std::string getDataFileNameByName (const std::string &relative_path) const
 
std::string getDataFilePath (const std::string &relative_path) const
 
const Parallel::Communicator & comm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static InputParameters validParams ()
 
static bool hasGenerateData (const InputParameters &params)
 
static void setHasGenerateData (InputParameters &params)
 

Public Attributes

const ConsoleStream _console
 

Static Public Attributes

static const std::string data_only_param
 
static constexpr auto SYSTEM
 
static constexpr auto NAME
 

Protected Member Functions

EChannelType getSubchannelType (unsigned int index) const
 returns the type of the subchannel given the index More...
 
Point rotatePoint (Point b, Real theta)
 
Point translatePoint (Point &b, Point &translation_vector)
 
Point getPinPosition (unsigned int i)
 returns the position of pin given pin index More...
 
std::vector< RealgetSubchannelPosition (unsigned int i)
 returns the position of subchannel given pin index More...
 
std::vector< unsigned intgetSubChannelPins (unsigned int i)
 returns the index of neighboring pins given subchannel index More...
 
virtual void generateData ()
 
T & copyMeshProperty (const std::string &target_data_name, const std::string &source_data_name, const std::string &source_mesh)
 
T & copyMeshProperty (const std::string &source_data_name, const std::string &source_mesh)
 
std::unique_ptr< MeshBase > & getMesh (const std::string &param_name, const bool allow_invalid=false)
 
std::vector< std::unique_ptr< MeshBase > *> getMeshes (const std::string &param_name)
 
std::unique_ptr< MeshBase > & getMeshByName (const MeshGeneratorName &mesh_generator_name)
 
std::vector< std::unique_ptr< MeshBase > *> getMeshesByName (const std::vector< MeshGeneratorName > &mesh_generator_names)
 
void declareMeshForSub (const std::string &param_name)
 
void declareMeshesForSub (const std::string &param_name)
 
void declareMeshForSubByName (const MeshGeneratorName &mesh_generator_name)
 
void declareMeshesForSubByName (const std::vector< MeshGeneratorName > &mesh_generator_names)
 
std::unique_ptr< MeshBase > buildMeshBaseObject (unsigned int dim=libMesh::invalid_uint)
 
std::unique_ptr< ReplicatedMesh > buildReplicatedMesh (unsigned int dim=libMesh::invalid_uint)
 
std::unique_ptr< DistributedMesh > buildDistributedMesh (unsigned int dim=libMesh::invalid_uint)
 
void addMeshSubgenerator (const std::string &type, const std::string &name, Ts... extra_input_parameters)
 
void addMeshSubgenerator (const std::string &type, const std::string &name, InputParameters params)
 
void declareNullMeshName (const MeshGeneratorName &name)
 
const T & getMeshProperty (const std::string &data_name, const std::string &prefix)
 
const T & getMeshProperty (const std::string &data_name)
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name, const std::string &prefix) const
 
bool hasMeshProperty (const std::string &data_name) const
 
bool hasMeshProperty (const std::string &data_name) const
 
std::string meshPropertyName (const std::string &data_name) const
 
T & declareMeshProperty (const std::string &data_name, Args &&... args)
 
T & declareMeshProperty (const std::string &data_name, const T &data_value)
 
T & declareMeshProperty (const std::string &data_name, Args &&... args)
 
T & declareMeshProperty (const std::string &data_name, const T &data_value)
 
T & setMeshProperty (const std::string &data_name, Args &&... args)
 
T & setMeshProperty (const std::string &data_name, const T &data_value)
 
T & setMeshProperty (const std::string &data_name, Args &&... args)
 
T & setMeshProperty (const std::string &data_name, const T &data_value)
 

Static Protected Member Functions

static std::string meshPropertyName (const std::string &data_name, const std::string &prefix)
 

Protected Attributes

const Real _unheated_length_entry
 unheated length of the fuel Pin at the entry of the assembly More...
 
const Real _heated_length
 heated length of the fuel Pin More...
 
const Real _unheated_length_exit
 unheated length of the fuel Pin at the exit of the assembly More...
 
std::vector< Real_z_grid
 axial location of nodes More...
 
const Real _pitch
 Distance between the neighbor fuel pins, pitch. More...
 
const Real _pin_diameter
 fuel Pin diameter More...
 
const unsigned int _n_rings
 Number of rings in the geometry. More...
 
const Real _flat_to_flat
 Half of gap between adjacent assemblies. More...
 
std::vector< EChannelType_subch_type
 Subchannel type. More...
 
std::vector< Point > _pin_position
 x,y coordinates of the fuel pins More...
 
std::vector< std::vector< Real > > _subchannel_position
 x,y coordinates of the subchannels More...
 
const unsigned int_block_id
 Subdomain ID used for the mesh block. More...
 
const unsigned int _n_cells
 Number of cells in the axial direction. More...
 
unsigned int _nrods
 Number of pins. More...
 
std::vector< std::vector< unsigned int > > _pins_in_rings
 fuel pins that are belonging to each ring More...
 
std::map< unsigned int, Real_orientation_map
 map inner and outer rings More...
 
unsigned int _n_channels
 number of subchannels More...
 
std::vector< std::vector< unsigned int > > _chan_to_pin_map
 stores the fuel pins belonging to each subchannel More...
 
const bool _verbose
 Flag to print out the detailed mesh assembly and coordinates. More...
 
MooseMesh *const _mesh
 
const bool & _enabled
 
MooseApp_app
 
const std::string _type
 
const std::string _name
 
const InputParameters_pars
 
Factory_factory
 
ActionFactory_action_factory
 
const Parallel::Communicator & _communicator
 

Detailed Description

Mesh generator that builds a 3D mesh representing triangular subchannels.

Definition at line 18 of file SCMDetailedTriSubChannelMeshGenerator.h.

Constructor & Destructor Documentation

◆ SCMDetailedTriSubChannelMeshGenerator()

SCMDetailedTriSubChannelMeshGenerator::SCMDetailedTriSubChannelMeshGenerator ( const InputParameters parameters)

Special case _n_rings == 1

Definition at line 43 of file SCMDetailedTriSubChannelMeshGenerator.C.

46  _unheated_length_entry(getParam<Real>("unheated_length_entry")),
47  _heated_length(getParam<Real>("heated_length")),
48  _unheated_length_exit(getParam<Real>("unheated_length_exit")),
49  _pitch(getParam<Real>("pitch")),
50  _pin_diameter(getParam<Real>("pin_diameter")),
51  _n_rings(getParam<unsigned int>("nrings")),
52  _flat_to_flat(getParam<Real>("flat_to_flat")),
53  _block_id(getParam<unsigned int>("block_id")),
54  _n_cells(getParam<unsigned int>("n_cells")),
55  _verbose(getParam<bool>("verbose_flag"))
56 {
58  Real dz = L / _n_cells;
59  for (unsigned int i = 0; i < _n_cells + 1; i++)
60  _z_grid.push_back(dz * i);
61 
62  // x coordinate for the first position
63  Real x0 = 0.0;
64  // y coordinate for the first position
65  Real y0 = 0.0;
66  // x coordinate for the second position
67  Real x1 = 0.0;
68  // y coordinate for the second position dummy variable
69  Real y1 = 0.0;
70  // dummy variable
71  Real a1 = 0.0;
72  // dummy variable
73  Real a2 = 0.0;
74  // average x coordinate
75  Real avg_coor_x = 0.0;
76  // average y coordinate
77  Real avg_coor_y = 0.0;
78  // distance between two points
79  Real dist = 0.0;
80  // distance between two points
81  Real dist0 = 0.0;
82  // the indicator used while setting _gap_to_chan_map array
83  std::vector<std::pair<unsigned int, unsigned int>> gap_fill;
85  _nrods = _pin_position.size();
86  // assign the pins to the corresponding rings
87  unsigned int k = 0; // initialize the fuel Pin counter index
88  _pins_in_rings.resize(_n_rings);
89  _pins_in_rings[0].push_back(k++);
90  for (unsigned int i = 1; i < _n_rings; i++)
91  for (unsigned int j = 0; j < i * 6; j++)
92  _pins_in_rings[i].push_back(k++);
93  // Given the number of pins and number of fuel Pin rings, the number of subchannels can be
94  // computed as follows:
95  unsigned int chancount = 0.0;
96  // Summing internal channels
97  for (unsigned int j = 0; j < _n_rings - 1; j++)
98  chancount += j * 6;
99  // Adding external channels to the total count
100  _n_channels = chancount + _nrods - 1 + (_n_rings - 1) * 6 + 6;
101 
102  // Utils for building the mesh
104  _subch_type.resize(_n_channels);
106 
107  for (unsigned int i = 0; i < _n_channels; i++)
108  {
109  _subchannel_position[i].reserve(3);
110  for (unsigned int j = 0; j < 3; j++)
111  {
112  _subchannel_position.at(i).push_back(0.0);
113  }
114  }
115 
116  // create the subchannels
117  k = 0; // initialize the subchannel counter index
118  for (unsigned int i = 1; i < _n_rings; i++)
119  {
120  // find the closest Pin at back ring
121  for (unsigned int j = 0; j < _pins_in_rings[i].size(); j++)
122  {
123  if (j == _pins_in_rings[i].size() - 1)
124  {
125  _chan_to_pin_map[k].push_back(_pins_in_rings[i][j]);
126  _chan_to_pin_map[k].push_back(_pins_in_rings[i][0]);
127  avg_coor_x =
128  0.5 * (_pin_position[_pins_in_rings[i][j]](0) + _pin_position[_pins_in_rings[i][0]](0));
129  avg_coor_y =
130  0.5 * (_pin_position[_pins_in_rings[i][j]](1) + _pin_position[_pins_in_rings[i][0]](1));
131  }
132  else
133  {
134  _chan_to_pin_map[k].push_back(_pins_in_rings[i][j]);
135  _chan_to_pin_map[k].push_back(_pins_in_rings[i][j + 1]);
136  avg_coor_x = 0.5 * (_pin_position[_pins_in_rings[i][j]](0) +
137  _pin_position[_pins_in_rings[i][j + 1]](0));
138  avg_coor_y = 0.5 * (_pin_position[_pins_in_rings[i][j]](1) +
139  _pin_position[_pins_in_rings[i][j + 1]](1));
140  }
141  dist0 = 1.0e+5;
142  _chan_to_pin_map[k].push_back(_pins_in_rings[i - 1][0]);
143  for (unsigned int l = 0; l < _pins_in_rings[i - 1].size(); l++)
144  {
145  dist = std::sqrt(pow(_pin_position[_pins_in_rings[i - 1][l]](0) - avg_coor_x, 2) +
146  pow(_pin_position[_pins_in_rings[i - 1][l]](1) - avg_coor_y, 2));
147 
148  if (dist < dist0)
149  {
150  _chan_to_pin_map[k][2] = _pins_in_rings[i - 1][l];
151  dist0 = dist;
152  }
153  }
155  _orientation_map.insert(std::make_pair(k, 0.0));
156  k = k + 1;
157  }
158 
159  // find the closest Pin at front ring
160  for (unsigned int j = 0; j < _pins_in_rings[i].size(); j++)
161  {
162  if (j == _pins_in_rings[i].size() - 1)
163  {
164  _chan_to_pin_map[k].push_back(_pins_in_rings[i][j]);
165  _chan_to_pin_map[k].push_back(_pins_in_rings[i][0]);
166  avg_coor_x =
167  0.5 * (_pin_position[_pins_in_rings[i][j]](0) + _pin_position[_pins_in_rings[i][0]](0));
168  avg_coor_y =
169  0.5 * (_pin_position[_pins_in_rings[i][j]](1) + _pin_position[_pins_in_rings[i][0]](1));
170  }
171  else
172  {
173  _chan_to_pin_map[k].push_back(_pins_in_rings[i][j]);
174  _chan_to_pin_map[k].push_back(_pins_in_rings[i][j + 1]);
175  avg_coor_x = 0.5 * (_pin_position[_pins_in_rings[i][j]](0) +
176  _pin_position[_pins_in_rings[i][j + 1]](0));
177  avg_coor_y = 0.5 * (_pin_position[_pins_in_rings[i][j]](1) +
178  _pin_position[_pins_in_rings[i][j + 1]](1));
179  }
180  // if the outermost ring, set the edge subchannels first... then the corner subchannels
181  if (i == _n_rings - 1)
182  {
183  // add edges
184  _subch_type[k] = EChannelType::EDGE; // an edge subchannel is created
185  k = k + 1;
186  if (j % i == 0)
187  {
188  // corner subchannel
189  _chan_to_pin_map[k].push_back(_pins_in_rings[i][j]);
191  k = k + 1;
192  }
193  // if not the outer most ring
194  }
195  else
196  {
197  dist0 = 1.0e+5;
198  _chan_to_pin_map[k].push_back(_pins_in_rings[i + 1][0]);
199  for (unsigned int l = 0; l < _pins_in_rings[i + 1].size(); l++)
200  {
201  dist = std::sqrt(pow(_pin_position[_pins_in_rings[i + 1][l]](0) - avg_coor_x, 2) +
202  pow(_pin_position[_pins_in_rings[i + 1][l]](1) - avg_coor_y, 2));
203  if (dist < dist0)
204  {
205  _chan_to_pin_map[k][2] = _pins_in_rings[i + 1][l];
206  dist0 = dist;
207  }
208  }
210  _orientation_map.insert(std::make_pair(k, libMesh::pi));
211  k = k + 1;
212  }
213  }
214  }
215 
216  for (auto & pin : _chan_to_pin_map)
217  pin.shrink_to_fit();
218 
219  // set the subchannel positions
220  Real _duct_to_pin_gap =
221  0.5 * (_flat_to_flat - (_n_rings - 1) * _pitch * std::sqrt(3.0) - _pin_diameter);
222  for (unsigned int i = 0; i < _n_channels; i++)
223  {
225  {
226  _subchannel_position[i][0] =
228  _pin_position[_chan_to_pin_map[i][2]](0)) /
229  3.0;
230  _subchannel_position[i][1] =
232  _pin_position[_chan_to_pin_map[i][2]](1)) /
233  3.0;
234  }
235  else if (_subch_type[i] == EChannelType::EDGE)
236  {
237  for (unsigned int j = 0; j < _n_channels; j++)
238  {
240  ((_chan_to_pin_map[i][0] == _chan_to_pin_map[j][0] &&
241  _chan_to_pin_map[i][1] == _chan_to_pin_map[j][1]) ||
242  (_chan_to_pin_map[i][0] == _chan_to_pin_map[j][1] &&
243  _chan_to_pin_map[i][1] == _chan_to_pin_map[j][0])))
244  {
245  x0 = _pin_position[_chan_to_pin_map[j][2]](0);
246  y0 = _pin_position[_chan_to_pin_map[j][2]](1);
247  }
248  else if (_subch_type[j] == EChannelType::CENTER &&
249  ((_chan_to_pin_map[i][0] == _chan_to_pin_map[j][0] &&
250  _chan_to_pin_map[i][1] == _chan_to_pin_map[j][2]) ||
251  (_chan_to_pin_map[i][0] == _chan_to_pin_map[j][2] &&
252  _chan_to_pin_map[i][1] == _chan_to_pin_map[j][0])))
253  {
254  x0 = _pin_position[_chan_to_pin_map[j][1]](0);
255  y0 = _pin_position[_chan_to_pin_map[j][1]](1);
256  }
257  else if (_subch_type[j] == EChannelType::CENTER &&
258  ((_chan_to_pin_map[i][0] == _chan_to_pin_map[j][1] &&
259  _chan_to_pin_map[i][1] == _chan_to_pin_map[j][2]) ||
260  (_chan_to_pin_map[i][0] == _chan_to_pin_map[j][2] &&
261  _chan_to_pin_map[i][1] == _chan_to_pin_map[j][1])))
262  {
263  x0 = _pin_position[_chan_to_pin_map[j][0]](0);
264  y0 = _pin_position[_chan_to_pin_map[j][0]](1);
265  }
266  x1 = 0.5 *
268  y1 = 0.5 *
270  a1 = _pin_diameter / 2.0 + _duct_to_pin_gap / 2.0;
271  a2 = std::sqrt((x1 - x0) * (x1 - x0) + (y1 - y0) * (y1 - y0)) + a1;
272  _subchannel_position[i][0] = (a2 * x1 - a1 * x0) / (a2 - a1);
273  _subchannel_position[i][1] = (a2 * y1 - a1 * y0) / (a2 - a1);
274  } // j
275  }
276  else if (_subch_type[i] == EChannelType::CORNER)
277  {
278  x0 = _pin_position[0](0);
279  y0 = _pin_position[0](1);
280  x1 = _pin_position[_chan_to_pin_map[i][0]](0);
281  y1 = _pin_position[_chan_to_pin_map[i][0]](1);
282  a1 = _pin_diameter / 2.0 + _duct_to_pin_gap / 2.0;
283  a2 = std::sqrt((x1 - x0) * (x1 - x0) + (y1 - y0) * (y1 - y0)) + a1;
284  _subchannel_position[i][0] = (a2 * x1 - a1 * x0) / (a2 - a1);
285  _subchannel_position[i][1] = (a2 * y1 - a1 * y0) / (a2 - a1);
286  }
287 
289  if (_n_rings == 1)
290  {
291  for (unsigned int i = 0; i < _n_channels; i++)
292  {
293  Real angle = (2 * i + 1) * libMesh::pi / 6.0;
295  _subchannel_position[i][0] = std::cos(angle) * _flat_to_flat / 2.0;
296  _subchannel_position[i][1] = std::sin(angle) * _flat_to_flat / 2.0;
297  }
298  }
299  }
300 }
const unsigned int _n_rings
Number of rings in the geometry.
const unsigned int & _block_id
Subdomain ID used for the mesh block.
std::vector< std::vector< Real > > _subchannel_position
x,y coordinates of the subchannels
const Real _pitch
Distance between the neighbor fuel pins, pitch.
MeshGenerator(const InputParameters &parameters)
std::vector< Real > _z_grid
axial location of nodes
const Real _unheated_length_exit
unheated length of the fuel Pin at the exit of the assembly
static void rodPositions(std::vector< Point > &positions, unsigned int nrings, Real pitch, Point center)
Calculates and stores the pin positions/centers for a hexagonal assembly containing the given number ...
ExpressionBuilder::EBTerm pow(const ExpressionBuilder::EBTerm &left, T exponent)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real _flat_to_flat
Half of gap between adjacent assemblies.
std::vector< Point > _pin_position
x,y coordinates of the fuel pins
std::vector< EChannelType > _subch_type
Subchannel type.
const InputParameters & parameters() const
std::vector< std::vector< unsigned int > > _pins_in_rings
fuel pins that are belonging to each ring
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const Real _unheated_length_entry
unheated length of the fuel Pin at the entry of the assembly
const unsigned int _n_cells
Number of cells in the axial direction.
std::vector< std::vector< unsigned int > > _chan_to_pin_map
stores the fuel pins belonging to each subchannel
const bool _verbose
Flag to print out the detailed mesh assembly and coordinates.
const Real _heated_length
heated length of the fuel Pin
static const std::string k
Definition: NS.h:130
const Real pi
std::map< unsigned int, Real > _orientation_map
map inner and outer rings

Member Function Documentation

◆ generate()

std::unique_ptr< MeshBase > SCMDetailedTriSubChannelMeshGenerator::generate ( )
overridevirtual

Implements MeshGenerator.

Definition at line 303 of file SCMDetailedTriSubChannelMeshGenerator.C.

304 {
305  auto mesh_base = buildMeshBaseObject();
306  BoundaryInfo & boundary_info = mesh_base->get_boundary_info();
307  mesh_base->set_spatial_dimension(3);
308  // Define the resolution (the number of points used to represent a circle).
309  // This must be divisible by 4.
310  const unsigned int theta_res_triangle = 24; // TODO: parameterize
311  const unsigned int theta_res_square = 24; // TODO: parameterize
312  // Compute the number of points needed to represent one sixth and quadrant of a circle.
313  const unsigned int points_per_sixth = theta_res_triangle / 6 + 1;
314  const unsigned int points_per_quadrant = theta_res_square / 4 + 1;
315 
316  // Compute the points needed to represent one axial cross-flow of a subchannel.
317  // For the center subchannel (sc) there is one center point plus the points from 3 side
318  // circles.
319  const unsigned int points_per_center = points_per_sixth * 3 + 1;
320  // For the corner sc there is one center point plus the points from 1 side circle plus 3
321  // corners
322  const unsigned int points_per_corner = points_per_sixth * 1 + 1 + 3;
323  // For the side sc there is one center point plus the points from 2 intersecting circles plus 2
324  // corners
325  const unsigned int points_per_side = points_per_quadrant * 2 + 1 + 2;
326 
327  if (_verbose)
328  {
329  _console << "Points per center: " << points_per_center << std::endl;
330  _console << "Points per side: " << points_per_side << std::endl;
331  _console << "Points per corner: " << points_per_corner << std::endl;
332  }
333 
334  // Compute the number of elements (Prism6) which combined base creates the sub-channel
335  // cross-section
336  const unsigned int elems_per_center = theta_res_triangle * 3 / 6 + 3; // TODO: check
337  const unsigned int elems_per_corner = theta_res_triangle / 6 + 4;
338  const unsigned int elems_per_side = 2 * theta_res_square / 4 + 4;
339  if (_verbose)
340  {
341  _console << "Elems per center: " << elems_per_center << std::endl;
342  _console << "Elems per side: " << elems_per_side << std::endl;
343  _console << "Elems per corner: " << elems_per_corner << std::endl;
344  }
345 
346  // specify number and type of sub-channel
347  unsigned int n_center, n_side, n_corner;
348  if (_n_rings == 1)
349  {
350  n_corner = 6;
351  n_side = 0;
352  n_center = _n_channels - n_side - n_corner;
353  }
354  else
355  {
356  n_corner = 6;
357  n_side = (_n_rings - 1) * 6;
358  n_center = _n_channels - n_side - n_corner;
359  }
360  if (_verbose)
361  {
362  _console << "Centers: " << n_center << std::endl;
363  _console << "Sides: " << n_side << std::endl;
364  _console << "Corners: " << n_corner << std::endl;
365  }
366 
367  // Compute the total number of points and elements.
368  const unsigned int points_per_level =
369  n_corner * points_per_corner + n_side * points_per_side + n_center * points_per_center;
370  const unsigned int elems_per_level =
371  n_corner * elems_per_corner + n_side * elems_per_side + n_center * elems_per_center;
372  if (_verbose)
373  {
374  _console << "Points per level: " << points_per_level << std::endl;
375  _console << "Elements per level: " << elems_per_level << std::endl;
376  }
377  const unsigned int n_points = points_per_level * (_n_cells + 1);
378  const unsigned int n_elems = elems_per_level * _n_cells;
379  if (_verbose)
380  {
381  _console << "Number of points: " << n_points << std::endl;
382  _console << "Number of elements: " << n_elems << std::endl;
383  }
384  mesh_base->reserve_nodes(n_points);
385  mesh_base->reserve_elem(n_elems);
386  // Build an array of points arranged in a circle on the xy-plane. (last and first node overlap)
387  // We build for both the square discretization in the edges and the triangular discretization
388  // within the mesh
389  const Real radius = _pin_diameter / 2.0;
390  std::array<Point, theta_res_square + 1> circle_points_square;
391  {
392  Real theta = 0;
393  for (unsigned int i = 0; i < theta_res_square + 1; i++)
394  {
395  circle_points_square[i](0) = radius * std::cos(theta);
396  circle_points_square[i](1) = radius * std::sin(theta);
397  theta += 2.0 * libMesh::pi / theta_res_square;
398  }
399  }
400  std::array<Point, theta_res_triangle + 1> circle_points_triangle;
401  {
402  Real theta = 0;
403  for (unsigned int i = 0; i < theta_res_triangle + 1; i++)
404  {
405  circle_points_triangle[i](0) = radius * std::cos(theta);
406  circle_points_triangle[i](1) = radius * std::sin(theta);
407  theta += 2.0 * libMesh::pi / theta_res_triangle;
408  }
409  }
410  // Define "quadrant center" reference points. These will be the centers of
411  // the 3 circles that represent the fuel pins. These centers are
412  // offset a little bit so that in the final mesh, there is a tiny gap between
413  // neighboring subchannel cells. That allows us to easily map a solution to
414  // this detailed mesh with a nearest-neighbor search.
415  const Real shrink_factor = 0.99999;
416  // Quadrants are used only for the side and corner subchannels
417  Real _duct_to_pin_gap =
418  0.5 * (_flat_to_flat - (_n_rings - 1) * _pitch * std::sqrt(3.0) - _pin_diameter);
419  std::array<Point, 2> quadrant_centers_sides;
420  quadrant_centers_sides[0] = Point(
421  -_pitch * 0.5 * shrink_factor, -(_duct_to_pin_gap + _pin_diameter) * 0.5 * shrink_factor, 0);
422  quadrant_centers_sides[1] = Point(
423  _pitch * 0.5 * shrink_factor, -(_duct_to_pin_gap + _pin_diameter) * 0.5 * shrink_factor, 0);
424  std::array<Point, 1> quadrant_centers_corner;
425  quadrant_centers_corner[0] =
426  Point(-(_duct_to_pin_gap + _pin_diameter) * 0.5 * std::sin(libMesh::pi / 6) * shrink_factor,
427  -(_duct_to_pin_gap + _pin_diameter) * 0.5 * std::cos(libMesh::pi / 6) * shrink_factor,
428  0);
429  // Triangles are used for all center subchannels
430  std::array<Point, 3> triangle_centers;
431  triangle_centers[0] = Point(0, _pitch * std::cos(libMesh::pi / 6) * 2 / 3 * shrink_factor, 0);
432  triangle_centers[1] = Point(-_pitch * 0.5 * shrink_factor,
433  -_pitch * std::cos(libMesh::pi / 6) * 1 / 3 * shrink_factor,
434  0);
435  triangle_centers[2] = Point(
436  _pitch * 0.5 * shrink_factor, -_pitch * std::cos(libMesh::pi / 6) * 1 / 3 * shrink_factor, 0);
437 
438  const unsigned int m_sixth = theta_res_triangle / 6;
439  const unsigned int m_quarter = theta_res_square / 4;
440  // Build an array of points that represent a cross section of a center subchannel
441  // cell. The points are ordered in this fashion:
442  // 3 1
443  // 2
444  // 0
445  // 4 5 8 9
446  // 6 7
447  std::array<Point, points_per_center> center_points;
448  {
449  unsigned int start;
450  for (unsigned int i = 0; i < 3; i++)
451  {
452  if (i == 0)
453  start = 5 * (m_sixth);
454  if (i == 1)
455  start = 1 * (m_sixth);
456  if (i == 2)
457  start = 3 * (m_sixth);
458  for (unsigned int ii = 0; ii < points_per_sixth; ii++)
459  {
460  auto c_pt = circle_points_triangle[start - ii];
461  center_points[i * points_per_sixth + ii + 1] = triangle_centers[i] + c_pt;
462  }
463  }
464  }
465 
466  // Build an array of points that represent a cross section of a top left corner subchannel
467  // cell. The points are ordered in this fashion (x direction towards point 4, y direction towards
468  // point 6):
469  // 5 6
470  //
471  // 0
472  // 4 2 3
473  // 1
474  std::array<Point, points_per_corner> corner_points;
475  {
476  for (unsigned int ii = 0; ii < points_per_sixth; ii++)
477  {
478  auto c_pt = circle_points_triangle[1 * m_quarter - ii];
479  corner_points[ii + 1] = quadrant_centers_corner[0] + c_pt;
480  }
481  Real side_short = (_duct_to_pin_gap + _pin_diameter) * 0.5;
482  Real side_long = (2.0 * _duct_to_pin_gap + _pin_diameter) * 0.5;
483  Real side_length = std::sqrt(std::pow(side_short, 2) + std::pow(side_long, 2) -
484  2 * side_short * side_long * std::cos(libMesh::pi / 6));
485  Real angle =
486  libMesh::pi - libMesh::pi / 3 -
487  std::acos((-std::pow(side_long, 2) + std::pow(side_short, 2) + std::pow(side_length, 2)) /
488  (2 * side_short * side_length));
489  corner_points[points_per_sixth + 1] = Point(side_length * std::cos(angle) * shrink_factor,
490  -side_length * std::sin(angle) * shrink_factor,
491  0);
492  corner_points[points_per_sixth + 2] =
493  Point(side_length * std::sin(libMesh::pi / 2 - angle - libMesh::pi / 6) * shrink_factor *
494  std::tan(libMesh::pi / 6),
495  side_length * std::sin(libMesh::pi / 2 - angle - libMesh::pi / 6) * shrink_factor,
496  0);
497  corner_points[points_per_sixth + 3] =
498  Point(-side_length * std::cos(libMesh::pi / 2 - angle - libMesh::pi / 6) * shrink_factor,
499  side_length * std::sin(libMesh::pi / 2 - angle - libMesh::pi / 6) * shrink_factor,
500  0);
501  }
502 
503  // Build an array of points that represent a cross-section of a top side subchannel
504  // cell. The points are ordered in this fashion:
505  // 8 7
506  //
507  // 0
508  // 1 2 5 6
509  // 3 4
510  std::array<Point, points_per_side> side_points;
511  {
512  for (unsigned int ii = 0; ii < points_per_quadrant; ii++)
513  {
514  auto c_pt = circle_points_square[m_quarter - ii];
515  side_points[ii + 1] = quadrant_centers_sides[0] + c_pt;
516  }
517  for (unsigned int ii = 0; ii < points_per_quadrant; ii++)
518  {
519  auto c_pt = circle_points_square[2 * m_quarter - ii];
520  side_points[points_per_quadrant + ii + 1] = quadrant_centers_sides[1] + c_pt;
521  }
522  side_points[2 * points_per_quadrant + 1] =
523  Point(_pitch * 0.5 * shrink_factor, 0.5 * _duct_to_pin_gap * shrink_factor, 0);
524  side_points[2 * points_per_quadrant + 2] =
525  Point(-_pitch * 0.5 * shrink_factor, 0.5 * _duct_to_pin_gap * shrink_factor, 0);
526  }
527 
528  int point_counter = 0;
529  unsigned int node_id = 0;
530  for (unsigned int i = 0; i < _n_channels; i++)
531  {
532  // Real offset_x = _flat_to_flat / 2.0;
533  // Real offset_y = _flat_to_flat / 2.0;
534 
536  {
537  for (auto z : _z_grid)
538  {
539  // Get height
540  Point z0{0, 0, z};
541 
542  // Get suchannel position and assign to point
543  auto loc_position = getSubchannelPosition(i);
544  Point p0{loc_position[0], loc_position[1], 0};
545 
546  // Determine orientation of current subchannel
547  auto subchannel_pins = getSubChannelPins(i);
548  Point subchannel_side =
549  getPinPosition(subchannel_pins[0]) + getPinPosition(subchannel_pins[1]);
550  Point base_center_orientation = {0, -1};
551 
552  // Get rotation angle for current subchannel
553  Real dot_prod = 0;
554  for (unsigned int lp = 0; lp < 2; lp++)
555  dot_prod += base_center_orientation(lp) * subchannel_side(lp);
556  auto theta =
557  std::acos(dot_prod / (base_center_orientation.norm() * subchannel_side.norm()));
558  if (subchannel_side(0) < 0)
559  theta = 2.0 * libMesh::pi - theta;
560 
561  // Real distance_side = subchannel_side.norm();
562  // Real distance_top = getPinPosition(subchannel_pins[2]).norm();
563  // if (distance_top > distance_side)
564  // theta += libMesh::pi * 0.0;
565 
566  theta += _orientation_map[i];
567 
568  theta = trunc((theta + (libMesh::pi / 6.0)) / (libMesh::pi / 3.0)) * libMesh::pi / 3.0;
569 
570  if (_verbose)
571  {
572  if (z == 0)
573  {
574  _console << "Subchannel Position: " << p0 << std::endl;
575  auto pins = getSubChannelPins(i);
576  for (auto r : pins)
577  _console << r << " ";
578  _console << std::endl;
579  _console << "Theta: " << theta / libMesh::pi * 180. << std::endl;
580  }
581  }
582 
583  // Assigning points for center channels
584  for (unsigned int i = 0; i < points_per_center; i++)
585  {
586  auto new_point = rotatePoint(center_points[i], theta) + p0 + z0;
587  if (_verbose)
588  {
589  if (z == 0)
590  _console << i << " - " << new_point << std::endl;
591  }
592  mesh_base->add_point(new_point, node_id++);
593  point_counter += 1;
594  }
595  }
596  }
597  else if (getSubchannelType(i) == EChannelType::EDGE)
598  {
599  for (auto z : _z_grid)
600  {
601  // Get height
602  Point z0{0, 0, z};
603 
604  // Get suchannel position and assign to point
605  auto loc_position = getSubchannelPosition(i);
606  Point p0{loc_position[0], loc_position[1], 0};
607 
608  // Determine orientation of current subchannel
609  auto subchannel_pins = getSubChannelPins(i);
610  Point subchannel_side =
611  getPinPosition(subchannel_pins[0]) + getPinPosition(subchannel_pins[1]);
612  Point base_center_orientation = {0, 1};
613 
614  // Get rotation angle for current subchannel
615  Real dot_prod = 0;
616  for (unsigned int lp = 0; lp < 2; lp++)
617  dot_prod += base_center_orientation(lp) * subchannel_side(lp);
618  auto theta =
619  std::acos(dot_prod / (base_center_orientation.norm() * subchannel_side.norm()));
620  if (subchannel_side(0) > 0)
621  theta = 2. * libMesh::pi - theta;
622  theta = trunc((theta + (libMesh::pi / 6.0)) / (libMesh::pi / 3.0)) * libMesh::pi / 3.0;
623 
624  if (_verbose)
625  {
626  if (z == 0)
627  {
628  _console << "Subchannel Position: " << p0 << std::endl;
629  auto pins = getSubChannelPins(i);
630  for (auto r : pins)
631  _console << r << " ";
632  _console << std::endl;
633  _console << "Theta: " << theta * 180 / libMesh::pi << std::endl;
634  }
635  }
636 
637  // Assigning points for center channels
638  for (unsigned int i = 0; i < points_per_side; i++)
639  {
640  auto new_point = rotatePoint(side_points[i], theta) + p0 + z0;
641  if (_verbose)
642  {
643  if (z == 0)
644  _console << i << " - " << new_point << std::endl;
645  }
646  mesh_base->add_point(new_point, node_id++);
647  point_counter += 1;
648  }
649  }
650  }
651  else // getSubchannelType(i) == EChannelType::CORNER
652  {
653  for (auto z : _z_grid)
654  {
655  // Get height
656  Point z0{0, 0, z};
657 
658  // Get suchannel position and assign to point
659  auto loc_position = getSubchannelPosition(i);
660  Point p0{loc_position[0], loc_position[1], 0};
661 
662  // Determine orientation of current subchannel
663  auto subchannel_pins = getSubChannelPins(i);
664  Point subchannel_side = getPinPosition(subchannel_pins[0]);
665  Point base_center_orientation = {1, 1};
666 
667  // Get rotation angle for current subchannel
668  Real dot_prod = 0;
669  for (unsigned int lp = 0; lp < 2; lp++)
670  dot_prod += base_center_orientation(lp) * subchannel_side(lp);
671  auto theta =
672  std::acos(dot_prod / (base_center_orientation.norm() * subchannel_side.norm()));
673  if (subchannel_side(0) > 0)
674  theta = 2. * libMesh::pi - theta;
675  theta = trunc((theta + (libMesh::pi / 6.0)) / (libMesh::pi / 3.0)) * libMesh::pi / 3.0;
676 
677  if (_verbose)
678  {
679  if (z == 0)
680  {
681  _console << "Subchannel Position: " << p0 << std::endl;
682  auto pins = getSubChannelPins(i);
683  for (auto r : pins)
684  _console << r << " ";
685  _console << std::endl;
686  _console << "Theta: " << theta * 180 / libMesh::pi << std::endl;
687  }
688  }
689 
690  // Assigning points for center channels
691  for (unsigned int i = 0; i < points_per_corner; i++)
692  {
693  auto new_point = rotatePoint(corner_points[i], theta) + p0 + z0;
694  if (_verbose)
695  {
696  if (z == 0)
697  _console << i << " - " << new_point << std::endl;
698  }
699  mesh_base->add_point(new_point, node_id++);
700  point_counter += 1;
701  }
702  }
703  }
704  } // i
705  if (_verbose)
706  _console << "Point counter: " << point_counter << std::endl;
707 
708  int element_counter = 0;
709  unsigned int elem_id = 0;
710  unsigned int number_of_corner = 0;
711  unsigned int number_of_side = 0;
712  unsigned int number_of_center = 0;
713  unsigned int elems_per_channel = 0;
714  unsigned int points_per_channel = 0;
715  for (unsigned int i = 0; i < _n_channels; i++)
716  {
717  auto subch_type = getSubchannelType(i);
718  if (subch_type == EChannelType::CORNER)
719  {
720  number_of_corner++;
721  elems_per_channel = elems_per_corner;
722  points_per_channel = points_per_corner;
723  if (_verbose)
724  _console << "Corner" << std::endl;
725  }
726  else if (subch_type == EChannelType::EDGE)
727  {
728  number_of_side++;
729  elems_per_channel = elems_per_side;
730  points_per_channel = points_per_side;
731  if (_verbose)
732  _console << "Edge" << std::endl;
733  }
734  else if (subch_type == EChannelType::CENTER)
735  {
736  number_of_center++;
737  elems_per_channel = elems_per_center;
738  points_per_channel = points_per_center;
739  if (_verbose)
740  _console << "Center" << std::endl;
741  }
742  for (unsigned int iz = 0; iz < _n_cells; iz++)
743  {
744  unsigned int elapsed_points = number_of_corner * points_per_corner * (_n_cells + 1) +
745  number_of_side * points_per_side * (_n_cells + 1) +
746  number_of_center * points_per_center * (_n_cells + 1) -
747  points_per_channel * (_n_cells + 1);
748  // index of the central node at base of cell
749  unsigned int indx1 = iz * points_per_channel + elapsed_points;
750  // index of the central node at top of cell
751  unsigned int indx2 = (iz + 1) * points_per_channel + elapsed_points;
752 
753  for (unsigned int i = 0; i < elems_per_channel; i++)
754  {
755  Elem * elem = new Prism6;
756  elem->subdomain_id() = _block_id;
757  elem->set_id(elem_id++);
758  elem = mesh_base->add_elem(elem);
759 
760  if (_verbose)
761  _console << "Node 0: " << *mesh_base->node_ptr(indx1) << std::endl;
762 
763  elem->set_node(0, mesh_base->node_ptr(indx1));
764  elem->set_node(1, mesh_base->node_ptr(indx1 + i + 1));
765  if (i != elems_per_channel - 1)
766  elem->set_node(2, mesh_base->node_ptr(indx1 + i + 2));
767  else
768  elem->set_node(2, mesh_base->node_ptr(indx1 + 1));
769 
770  elem->set_node(3, mesh_base->node_ptr(indx2));
771  elem->set_node(4, mesh_base->node_ptr(indx2 + i + 1));
772  if (i != elems_per_channel - 1)
773  elem->set_node(5, mesh_base->node_ptr(indx2 + i + 2));
774  else
775  elem->set_node(5, mesh_base->node_ptr(indx2 + 1));
776 
777  if (iz == 0)
778  boundary_info.add_side(elem, 0, 0);
779  if (iz == _n_cells - 1)
780  boundary_info.add_side(elem, 4, 1);
781 
782  element_counter += 1;
783  }
784  }
785  }
786  if (_verbose)
787  _console << "Element counter: " << element_counter << std::endl;
788  boundary_info.sideset_name(0) = "inlet";
789  boundary_info.sideset_name(1) = "outlet";
790  mesh_base->subdomain_name(_block_id) = name();
791  if (_verbose)
792  _console << "Mesh assembly done" << std::endl;
793  mesh_base->prepare_for_use();
794 
795  return mesh_base;
796 }
const unsigned int _n_rings
Number of rings in the geometry.
const unsigned int & _block_id
Subdomain ID used for the mesh block.
const Real radius
const Real _pitch
Distance between the neighbor fuel pins, pitch.
Point getPinPosition(unsigned int i)
returns the position of pin given pin index
EChannelType getSubchannelType(unsigned int index) const
returns the type of the subchannel given the index
std::vector< Real > _z_grid
axial location of nodes
virtual const std::string & name() const
std::vector< unsigned int > getSubChannelPins(unsigned int i)
returns the index of neighboring pins given subchannel index
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const Real _flat_to_flat
Half of gap between adjacent assemblies.
std::vector< Real > getSubchannelPosition(unsigned int i)
returns the position of subchannel given pin index
const ConsoleStream _console
std::unique_ptr< MeshBase > buildMeshBaseObject(unsigned int dim=libMesh::invalid_uint)
const unsigned int _n_cells
Number of cells in the axial direction.
const bool _verbose
Flag to print out the detailed mesh assembly and coordinates.
MooseUnits pow(const MooseUnits &, int)
const Real pi
std::map< unsigned int, Real > _orientation_map
map inner and outer rings

◆ getPinPosition()

Point SCMDetailedTriSubChannelMeshGenerator::getPinPosition ( unsigned int  i)
inlineprotected

returns the position of pin given pin index

Definition at line 30 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate().

30 { return _pin_position[i]; }
std::vector< Point > _pin_position
x,y coordinates of the fuel pins

◆ getSubChannelPins()

std::vector<unsigned int> SCMDetailedTriSubChannelMeshGenerator::getSubChannelPins ( unsigned int  i)
inlineprotected

returns the index of neighboring pins given subchannel index

Definition at line 34 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate().

34 { return _chan_to_pin_map[i]; }
std::vector< std::vector< unsigned int > > _chan_to_pin_map
stores the fuel pins belonging to each subchannel

◆ getSubchannelPosition()

std::vector<Real> SCMDetailedTriSubChannelMeshGenerator::getSubchannelPosition ( unsigned int  i)
inlineprotected

returns the position of subchannel given pin index

Definition at line 32 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate().

32 { return _subchannel_position[i]; }
std::vector< std::vector< Real > > _subchannel_position
x,y coordinates of the subchannels

◆ getSubchannelType()

EChannelType SCMDetailedTriSubChannelMeshGenerator::getSubchannelType ( unsigned int  index) const
inlineprotected

returns the type of the subchannel given the index

Definition at line 26 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate().

26 { return _subch_type[index]; }
std::vector< EChannelType > _subch_type
Subchannel type.

◆ rotatePoint()

Point SCMDetailedTriSubChannelMeshGenerator::rotatePoint ( Point  b,
Real  theta 
)
protected

Definition at line 799 of file SCMDetailedTriSubChannelMeshGenerator.C.

Referenced by generate().

800 {
801 
802  // Building rotation matrix
803  std::vector<std::vector<Real>> A;
804  A.resize(3);
805  for (std::vector<Real> a : A)
806  {
807  a.resize(3);
808  }
809 
810  A[0] = {std::cos(theta), -std::sin(theta), 0.0};
811  A[1] = {std::sin(theta), std::cos(theta), 0.0};
812  A[2] = {0.0, 0.0, 1.0};
813 
814  // Rotating vector
815  Point rotated_vector = Point(0.0, 0.0, 0.0);
816  for (unsigned int i = 0; i < 3; i++)
817  {
818  for (unsigned int j = 0; j < 3; j++)
819  {
820  rotated_vector(i) += A[i][j] * b(j);
821  }
822  }
823 
824  return rotated_vector;
825 }
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")

◆ translatePoint()

Point SCMDetailedTriSubChannelMeshGenerator::translatePoint ( Point &  b,
Point &  translation_vector 
)
protected

Definition at line 828 of file SCMDetailedTriSubChannelMeshGenerator.C.

829 {
830  // Translating point
831  Point translated_vector = Point(0.0, 0.0, 0.0);
832  for (unsigned int i = 0; i < 3; i++)
833  {
834  translated_vector(i) = b(i) + translation_vector(i);
835  }
836 
837  return translated_vector;
838 }

◆ validParams()

InputParameters SCMDetailedTriSubChannelMeshGenerator::validParams ( )
static

Definition at line 24 of file SCMDetailedTriSubChannelMeshGenerator.C.

25 {
27  params.addClassDescription(
28  "Creates a detailed mesh of subchannels in a triangular lattice arrangement");
29  params.addRequiredParam<Real>("pitch", "Pitch [m]");
30  params.addRequiredParam<Real>("pin_diameter", "Rod diameter [m]");
31  params.addParam<Real>("unheated_length_entry", 0.0, "Unheated length at entry [m]");
32  params.addRequiredParam<Real>("heated_length", "Heated length [m]");
33  params.addParam<Real>("unheated_length_exit", 0.0, "Unheated length at exit [m]");
34  params.addRequiredParam<unsigned int>("nrings", "Number of fuel Pin rings per assembly [-]");
35  params.addRequiredParam<Real>("flat_to_flat",
36  "Flat to flat distance for the hexagonal assembly [m]");
37  params.addParam<unsigned int>("block_id", 0, "Block ID used for the mesh subdomain.");
38  params.addRequiredParam<unsigned int>("n_cells", "The number of cells in the axial direction");
39  params.addParam<bool>("verbose_flag", false, "Flag to print out the mesh coordinates.");
40  return params;
41 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void addRequiredParam(const std::string &name, const std::string &doc_string)
static InputParameters validParams()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _block_id

const unsigned int& SCMDetailedTriSubChannelMeshGenerator::_block_id
protected

Subdomain ID used for the mesh block.

Definition at line 59 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate().

◆ _chan_to_pin_map

std::vector<std::vector<unsigned int> > SCMDetailedTriSubChannelMeshGenerator::_chan_to_pin_map
protected

stores the fuel pins belonging to each subchannel

Definition at line 71 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by getSubChannelPins(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _flat_to_flat

const Real SCMDetailedTriSubChannelMeshGenerator::_flat_to_flat
protected

Half of gap between adjacent assemblies.

Definition at line 51 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _heated_length

const Real SCMDetailedTriSubChannelMeshGenerator::_heated_length
protected

heated length of the fuel Pin

Definition at line 39 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by SCMDetailedTriSubChannelMeshGenerator().

◆ _n_cells

const unsigned int SCMDetailedTriSubChannelMeshGenerator::_n_cells
protected

Number of cells in the axial direction.

Definition at line 61 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _n_channels

unsigned int SCMDetailedTriSubChannelMeshGenerator::_n_channels
protected

number of subchannels

Definition at line 69 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _n_rings

const unsigned int SCMDetailedTriSubChannelMeshGenerator::_n_rings
protected

Number of rings in the geometry.

Definition at line 49 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _nrods

unsigned int SCMDetailedTriSubChannelMeshGenerator::_nrods
protected

Number of pins.

Definition at line 63 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by SCMDetailedTriSubChannelMeshGenerator().

◆ _orientation_map

std::map<unsigned int, Real> SCMDetailedTriSubChannelMeshGenerator::_orientation_map
protected

map inner and outer rings

Definition at line 67 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _pin_diameter

const Real SCMDetailedTriSubChannelMeshGenerator::_pin_diameter
protected

fuel Pin diameter

Definition at line 47 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _pin_position

std::vector<Point> SCMDetailedTriSubChannelMeshGenerator::_pin_position
protected

x,y coordinates of the fuel pins

Definition at line 55 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by getPinPosition(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _pins_in_rings

std::vector<std::vector<unsigned int> > SCMDetailedTriSubChannelMeshGenerator::_pins_in_rings
protected

fuel pins that are belonging to each ring

Definition at line 65 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by SCMDetailedTriSubChannelMeshGenerator().

◆ _pitch

const Real SCMDetailedTriSubChannelMeshGenerator::_pitch
protected

Distance between the neighbor fuel pins, pitch.

Definition at line 45 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _subch_type

std::vector<EChannelType> SCMDetailedTriSubChannelMeshGenerator::_subch_type
protected

Subchannel type.

Definition at line 53 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by getSubchannelType(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _subchannel_position

std::vector<std::vector<Real> > SCMDetailedTriSubChannelMeshGenerator::_subchannel_position
protected

x,y coordinates of the subchannels

Definition at line 57 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by getSubchannelPosition(), and SCMDetailedTriSubChannelMeshGenerator().

◆ _unheated_length_entry

const Real SCMDetailedTriSubChannelMeshGenerator::_unheated_length_entry
protected

unheated length of the fuel Pin at the entry of the assembly

Definition at line 37 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by SCMDetailedTriSubChannelMeshGenerator().

◆ _unheated_length_exit

const Real SCMDetailedTriSubChannelMeshGenerator::_unheated_length_exit
protected

unheated length of the fuel Pin at the exit of the assembly

Definition at line 41 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by SCMDetailedTriSubChannelMeshGenerator().

◆ _verbose

const bool SCMDetailedTriSubChannelMeshGenerator::_verbose
protected

Flag to print out the detailed mesh assembly and coordinates.

Definition at line 73 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate().

◆ _z_grid

std::vector<Real> SCMDetailedTriSubChannelMeshGenerator::_z_grid
protected

axial location of nodes

Definition at line 43 of file SCMDetailedTriSubChannelMeshGenerator.h.

Referenced by generate(), and SCMDetailedTriSubChannelMeshGenerator().


The documentation for this class was generated from the following files: