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

Defines a discretized line segment in 3D space. More...

#include <DiscreteLineSegmentInterface.h>

Inheritance diagram for DiscreteLineSegmentInterface:
[legend]

Public Member Functions

 DiscreteLineSegmentInterface (const MooseObject *moose_object)
 
virtual Point getPosition () const
 
virtual RealVectorValue getDirection () const
 
virtual Real getRotation () const
 
virtual Real getNumElems () const
 
virtual Real getLength () const
 
Real computeAxialCoordinate (const Point &p) const
 
Real computeRadialCoordinate (const Point &p) const
 
unsigned int getAxialSectionIndex (const Point &p) const
 
unsigned int getAxialElementIndex (const Point &p_center) const
 
Point computeRealPointFromReferencePoint (const Point &p) const
 Computes point in 3-D space from a point in reference space. More...
 
Point computeReferencePointFromRealPoint (const Point &p) const
 Computes point in reference space from a point in 3-D space. More...
 
MooseEnum getAlignmentAxis () const
 Gets an axis MooseEnum for the axis the component is aligned with. More...
 
std::vector< RealgetElementBoundaryCoordinates () const
 Gets the element boundary coordinates for the aligned axis. More...
 

Static Public Member Functions

static InputParameters validParams ()
 
static RealTensorValue computeDirectionTransformationTensor (const RealVectorValue &dir)
 Computes the direction transformation tensor. More...
 
static RealTensorValue computeXRotationTransformationTensor (const Real &rotation)
 Computes the rotation transformation tensor. More...
 
static Point computeRealPointFromReferencePoint (const Point &p, const RealVectorValue &position, const RealTensorValue &R, const RealTensorValue &Rx)
 Computes point in 3-D space from a point in reference space. More...
 
static Point computeReferencePointFromRealPoint (const Point &p, const RealVectorValue &position, const RealTensorValue &R_inv, const RealTensorValue &Rx_inv)
 Computes point in reference space from a point in 3-D space. More...
 
static MooseEnum getAlignmentAxis (const RealVectorValue &dir)
 Gets an axis MooseEnum for the axis the component is aligned with. More...
 
static std::vector< RealgetElementBoundaryCoordinates (const RealVectorValue &position, const RealVectorValue &orientation, const Real &rotation, const std::vector< Real > &lengths, const std::vector< unsigned int > &n_elems)
 Gets the element boundary coordinates for the aligned axis. More...
 

Protected Attributes

const Point & _position
 Start position of axis in 3-D space. More...
 
const RealVectorValue_dir_unnormalized
 Unnormalized direction of axis from start position to end position. More...
 
const RealVectorValue _dir
 Normalized direction of axis from start position to end position. More...
 
const Real_rotation
 Angle of rotation about the x-axis. More...
 
std::vector< Real_lengths
 Length of each axial section. More...
 
Real _length
 Total axial length. More...
 
const std::vector< unsigned int > & _n_elems
 Number of elements in each axial section. More...
 
const unsigned int _n_elem
 Total number of axial elements. More...
 
const unsigned int _n_sections
 Number of axial sections. More...
 
std::vector< Real_section_end
 Axial coordinate of the end of each axial section using the line 'position' as the origin. More...
 
std::vector< Real_x_centers
 Center axial coordinate of each axial element. More...
 
const RealTensorValue _R
 Direction transformation tensor. More...
 
const RealTensorValue _Rx
 Rotational transformation tensor about x-axis. More...
 
const RealTensorValue _R_inv
 Inverse direction transformation tensor. More...
 
const RealTensorValue _Rx_inv
 Inverse rotational transformation tensor about x-axis. More...
 
const std::string _moose_object_name_dlsi
 Name of the MOOSE object. More...
 

Static Private Member Functions

static RealVectorValue initializeDirectionVector (const RealVectorValue &dir_unnormalized)
 Computes a normalized direction vector or reports an error if the zero vector is provided. More...
 

Detailed Description

Defines a discretized line segment in 3D space.

Definition at line 22 of file DiscreteLineSegmentInterface.h.

Constructor & Destructor Documentation

◆ DiscreteLineSegmentInterface()

DiscreteLineSegmentInterface::DiscreteLineSegmentInterface ( const MooseObject moose_object)

Definition at line 33 of file DiscreteLineSegmentInterface.C.

34  : _position(moose_object->parameters().get<Point>("position")),
35  _dir_unnormalized(moose_object->parameters().get<RealVectorValue>("orientation")),
37  _rotation(moose_object->parameters().get<Real>("rotation")),
38  _lengths(moose_object->parameters().get<std::vector<Real>>("length")),
39  _length(std::accumulate(_lengths.begin(), _lengths.end(), 0.0)),
40  _n_elems(moose_object->parameters().get<std::vector<unsigned int>>("n_elems")),
41  _n_elem(std::accumulate(_n_elems.begin(), _n_elems.end(), 0)),
42  _n_sections(_lengths.size()),
47  _R_inv(_R.inverse()),
48  _Rx_inv(_Rx.inverse()),
49  _moose_object_name_dlsi(moose_object->name())
50 {
51  std::partial_sum(_lengths.begin(), _lengths.end(), _section_end.begin());
52 
53  if (_lengths.size() != _n_elems.size())
54  mooseError("The parameters 'length' and 'n_elems' must have the same number of entries.");
55 
56  // Compute the axial coordinates of the centers of each element
57  unsigned int k_section_begin = 0;
58  Real x_begin = 0.0;
59  for (unsigned int j = 0; j < _n_sections; j++)
60  {
61  const Real dx = _lengths[j] / _n_elems[j];
62  for (unsigned int i = 0; i < _n_elems[j]; i++)
63  {
64  const unsigned int k = k_section_begin + i;
65  _x_centers[k] = x_begin + 0.5 * dx;
66  x_begin += dx;
67  }
68  k_section_begin += _n_elems[j];
69  }
70 }
void mooseError(Args &&... args)
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
static RealTensorValue computeXRotationTransformationTensor(const Real &rotation)
Computes the rotation transformation tensor.
const RealTensorValue _Rx_inv
Inverse rotational transformation tensor about x-axis.
std::vector< Real > _lengths
Length of each axial section.
static RealTensorValue computeDirectionTransformationTensor(const RealVectorValue &dir)
Computes the direction transformation tensor.
const Real & _rotation
Angle of rotation about the x-axis.
const RealVectorValue & _dir_unnormalized
Unnormalized direction of axis from start position to end position.
virtual const std::string & name() const
const RealVectorValue _dir
Normalized direction of axis from start position to end position.
std::vector< Real > _section_end
Axial coordinate of the end of each axial section using the line &#39;position&#39; as the origin...
const unsigned int _n_sections
Number of axial sections.
const Point & _position
Start position of axis in 3-D space.
const std::vector< unsigned int > & _n_elems
Number of elements in each axial section.
const unsigned int _n_elem
Total number of axial elements.
const RealTensorValue _R_inv
Inverse direction transformation tensor.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const RealTensorValue _Rx
Rotational transformation tensor about x-axis.
const RealTensorValue _R
Direction transformation tensor.
const InputParameters & parameters() const
const std::string _moose_object_name_dlsi
Name of the MOOSE object.
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
static RealVectorValue initializeDirectionVector(const RealVectorValue &dir_unnormalized)
Computes a normalized direction vector or reports an error if the zero vector is provided.
static const std::string k
Definition: NS.h:130
std::vector< Real > _x_centers
Center axial coordinate of each axial element.

Member Function Documentation

◆ computeAxialCoordinate()

Real DiscreteLineSegmentInterface::computeAxialCoordinate ( const Point &  p) const

Definition at line 107 of file DiscreteLineSegmentInterface.C.

Referenced by DiscreteLineSegmentInterfaceTestAux::computeValue(), getAxialElementIndex(), and getAxialSectionIndex().

108 {
109  const Real ax_coord = _dir * (p - _position);
110  if (MooseUtils::absoluteFuzzyLessThan(ax_coord, 0.0) ||
113  ": The point ",
114  p,
115  " has an invalid axial coordinate (",
116  ax_coord,
117  "). Valid axial coordinates are in the range (0,",
118  _length,
119  ").");
120  else
121  return ax_coord;
122 }
void mooseError(Args &&... args)
const RealVectorValue _dir
Normalized direction of axis from start position to end position.
bool absoluteFuzzyLessThan(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
const Point & _position
Start position of axis in 3-D space.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const std::string _moose_object_name_dlsi
Name of the MOOSE object.
bool absoluteFuzzyGreaterThan(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)

◆ computeDirectionTransformationTensor()

RealTensorValue DiscreteLineSegmentInterface::computeDirectionTransformationTensor ( const RealVectorValue dir)
static

Computes the direction transformation tensor.

Parameters
[in]dirDirection vector

Definition at line 82 of file DiscreteLineSegmentInterface.C.

Referenced by getElementBoundaryCoordinates().

83 {
84  if (MooseUtils::absoluteFuzzyEqual(dir.norm(), 0.0))
85  mooseError("The direction vector must not be the zero vector.");
86 
87  const auto dir_normalized = dir / dir.norm();
88  const Real theta = acos(dir_normalized(2));
89  const Real aphi = atan2(dir_normalized(1), dir_normalized(0));
90  return RealTensorValue(
91  RealVectorValue(cos(aphi) * sin(theta), -sin(aphi), -cos(aphi) * cos(theta)),
92  RealVectorValue(sin(aphi) * sin(theta), cos(aphi), -sin(aphi) * cos(theta)),
93  RealVectorValue(cos(theta), 0.0, sin(theta)));
94 }
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sin(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tan
auto norm() const -> decltype(std::norm(Real()))
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
void mooseError(Args &&... args)
TensorValue< Real > RealTensorValue
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template cos(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(cos
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ computeRadialCoordinate()

Real DiscreteLineSegmentInterface::computeRadialCoordinate ( const Point &  p) const

Definition at line 125 of file DiscreteLineSegmentInterface.C.

Referenced by DiscreteLineSegmentInterfaceTestAux::computeValue().

126 {
127  const RealVectorValue v = (p - _position);
128  return v.cross(_dir).norm();
129 }
const RealVectorValue _dir
Normalized direction of axis from start position to end position.
const Point & _position
Start position of axis in 3-D space.
static const std::string v
Definition: NS.h:84

◆ computeRealPointFromReferencePoint() [1/2]

Point DiscreteLineSegmentInterface::computeRealPointFromReferencePoint ( const Point &  p) const

Computes point in 3-D space from a point in reference space.

Parameters
[in]pPoint in reference space

Definition at line 163 of file DiscreteLineSegmentInterface.C.

Referenced by getElementBoundaryCoordinates(), and GeneratedMeshComponent::setupMesh().

164 {
166 }
const Point & _position
Start position of axis in 3-D space.
const RealTensorValue _Rx
Rotational transformation tensor about x-axis.
const RealTensorValue _R
Direction transformation tensor.
Point computeRealPointFromReferencePoint(const Point &p) const
Computes point in 3-D space from a point in reference space.

◆ computeRealPointFromReferencePoint() [2/2]

Point DiscreteLineSegmentInterface::computeRealPointFromReferencePoint ( const Point &  p,
const RealVectorValue position,
const RealTensorValue R,
const RealTensorValue Rx 
)
static

Computes point in 3-D space from a point in reference space.

Parameters
[in]pPoint in reference space

Definition at line 154 of file DiscreteLineSegmentInterface.C.

158 {
159  return R * (Rx * p) + position;
160 }
static const std::string R
Definition: NS.h:162

◆ computeReferencePointFromRealPoint() [1/2]

Point DiscreteLineSegmentInterface::computeReferencePointFromRealPoint ( const Point &  p) const

Computes point in reference space from a point in 3-D space.

Parameters
[in]pPoint in 3-D space

Definition at line 178 of file DiscreteLineSegmentInterface.C.

179 {
181 }
const RealTensorValue _Rx_inv
Inverse rotational transformation tensor about x-axis.
const Point & _position
Start position of axis in 3-D space.
Point computeReferencePointFromRealPoint(const Point &p) const
Computes point in reference space from a point in 3-D space.
const RealTensorValue _R_inv
Inverse direction transformation tensor.

◆ computeReferencePointFromRealPoint() [2/2]

Point DiscreteLineSegmentInterface::computeReferencePointFromRealPoint ( const Point &  p,
const RealVectorValue position,
const RealTensorValue R_inv,
const RealTensorValue Rx_inv 
)
static

Computes point in reference space from a point in 3-D space.

Parameters
[in]pPoint in 3-D space

Definition at line 169 of file DiscreteLineSegmentInterface.C.

173 {
174  return Rx_inv * R_inv * (p - position);
175 }

◆ computeXRotationTransformationTensor()

RealTensorValue DiscreteLineSegmentInterface::computeXRotationTransformationTensor ( const Real rotation)
static

Computes the rotation transformation tensor.

Parameters
[in]rotationRotation about the x-axis, in degrees

Definition at line 97 of file DiscreteLineSegmentInterface.C.

Referenced by getElementBoundaryCoordinates().

98 {
99  const Real rotation_rad = M_PI * rotation / 180.;
100  const RealVectorValue Rx_x(1., 0., 0.);
101  const RealVectorValue Rx_y(0., cos(rotation_rad), -sin(rotation_rad));
102  const RealVectorValue Rx_z(0., sin(rotation_rad), cos(rotation_rad));
103  return RealTensorValue(Rx_x, Rx_y, Rx_z);
104 }
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template * sin(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(tan
TensorValue< Real > RealTensorValue
CTSub CT_OPERATOR_BINARY CTMul CTCompareLess CTCompareGreater CTCompareEqual _arg template cos(_arg) *_arg.template D< dtag >()) CT_SIMPLE_UNARY_FUNCTION(cos
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ getAlignmentAxis() [1/2]

MooseEnum DiscreteLineSegmentInterface::getAlignmentAxis ( ) const

Gets an axis MooseEnum for the axis the component is aligned with.

If the axis does not align with the x, y, or z axis, then an invalid MooseEnum is returned.

Definition at line 198 of file DiscreteLineSegmentInterface.C.

Referenced by CoupledHeatTransferAction::addUserObjects(), CoupledHeatTransferAction::CoupledHeatTransferAction(), and getElementBoundaryCoordinates().

199 {
200  return getAlignmentAxis(_dir);
201 }
const RealVectorValue _dir
Normalized direction of axis from start position to end position.
MooseEnum getAlignmentAxis() const
Gets an axis MooseEnum for the axis the component is aligned with.

◆ getAlignmentAxis() [2/2]

MooseEnum DiscreteLineSegmentInterface::getAlignmentAxis ( const RealVectorValue dir)
static

Gets an axis MooseEnum for the axis the component is aligned with.

If the axis does not align with the x, y, or z axis, then an invalid MooseEnum is returned.

Definition at line 184 of file DiscreteLineSegmentInterface.C.

185 {
186  MooseEnum axis("x y z");
187  if (THM::areParallelVectors(dir, RealVectorValue(1, 0, 0)))
188  axis = "x";
189  else if (THM::areParallelVectors(dir, RealVectorValue(0, 1, 0)))
190  axis = "y";
191  else if (THM::areParallelVectors(dir, RealVectorValue(0, 0, 1)))
192  axis = "z";
193 
194  return axis;
195 }
static const std::string axis
Definition: NS.h:27
bool areParallelVectors(const RealVectorValue &a, const RealVectorValue &b, const Real &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
Tests if two real-valued vectors are parallel within some absolute tolerance.
Definition: Numerics.C:26

◆ getAxialElementIndex()

unsigned int DiscreteLineSegmentInterface::getAxialElementIndex ( const Point &  p_center) const

Definition at line 143 of file DiscreteLineSegmentInterface.C.

Referenced by DiscreteLineSegmentInterfaceTestAux::computeValue().

144 {
145  const Real axial_coordinate = computeAxialCoordinate(p_center);
146  for (unsigned int i = 0; i < _n_elem; ++i)
147  if (MooseUtils::absoluteFuzzyEqual(axial_coordinate, _x_centers[i]))
148  return i;
149 
150  mooseError("No axial element index was found.");
151 }
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
void mooseError(Args &&... args)
Real computeAxialCoordinate(const Point &p) const
const unsigned int _n_elem
Total number of axial elements.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< Real > _x_centers
Center axial coordinate of each axial element.

◆ getAxialSectionIndex()

unsigned int DiscreteLineSegmentInterface::getAxialSectionIndex ( const Point &  p) const

Definition at line 132 of file DiscreteLineSegmentInterface.C.

Referenced by DiscreteLineSegmentInterfaceTestAux::computeValue().

133 {
134  const Real axial_coordinate = computeAxialCoordinate(p);
135  for (unsigned int i = 0; i < _n_sections; ++i)
136  if (MooseUtils::absoluteFuzzyLessEqual(axial_coordinate, _section_end[i]))
137  return i;
138 
139  mooseError("No axial section index was found.");
140 }
void mooseError(Args &&... args)
Real computeAxialCoordinate(const Point &p) const
std::vector< Real > _section_end
Axial coordinate of the end of each axial section using the line &#39;position&#39; as the origin...
const unsigned int _n_sections
Number of axial sections.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool absoluteFuzzyLessEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)

◆ getDirection()

virtual RealVectorValue DiscreteLineSegmentInterface::getDirection ( ) const
inlinevirtual

◆ getElementBoundaryCoordinates() [1/2]

std::vector< Real > DiscreteLineSegmentInterface::getElementBoundaryCoordinates ( ) const

Gets the element boundary coordinates for the aligned axis.

Definition at line 254 of file DiscreteLineSegmentInterface.C.

Referenced by CoupledHeatTransferAction::addUserObjects().

255 {
257 }
std::vector< Real > _lengths
Length of each axial section.
const Real & _rotation
Angle of rotation about the x-axis.
const RealVectorValue _dir
Normalized direction of axis from start position to end position.
const Point & _position
Start position of axis in 3-D space.
const std::vector< unsigned int > & _n_elems
Number of elements in each axial section.
std::vector< Real > getElementBoundaryCoordinates() const
Gets the element boundary coordinates for the aligned axis.

◆ getElementBoundaryCoordinates() [2/2]

std::vector< Real > DiscreteLineSegmentInterface::getElementBoundaryCoordinates ( const RealVectorValue position,
const RealVectorValue orientation,
const Real rotation,
const std::vector< Real > &  lengths,
const std::vector< unsigned int > &  n_elems 
)
static

Gets the element boundary coordinates for the aligned axis.

Parameters
[in]positionStart position of axis in 3-D space
[in]orientationDirection of axis from start position to end position
[in]rotationAngle of rotation about the x-axis, in degrees
[in]lengthsLength of each axial section
[in]n_elemsNumber of elements in each axial section

Definition at line 204 of file DiscreteLineSegmentInterface.C.

210 {
211  const auto axis = getAlignmentAxis(orientation);
212 
213  unsigned int d;
214  if (axis == "x")
215  d = 0;
216  else if (axis == "y")
217  d = 1;
218  else if (axis == "z")
219  d = 2;
220  else
221  mooseError("Invalid axis.");
222 
223  const auto R = computeDirectionTransformationTensor(orientation);
224  const auto Rx = computeXRotationTransformationTensor(rotation);
225 
226  const unsigned int n_elems_total = std::accumulate(n_elems.begin(), n_elems.end(), 0);
227  const unsigned int n_sections = lengths.size();
228 
229  unsigned int i_section_start = 0;
230  Real section_start_ref_position = 0.0;
231  std::vector<Real> element_boundary_coordinates(n_elems_total + 1);
232  element_boundary_coordinates[0] =
233  computeRealPointFromReferencePoint(Point(0, 0, 0), position, R, Rx)(d);
234  for (unsigned int j = 0; j < n_sections; ++j)
235  {
236  const Real section_dx = lengths[j] / n_elems[j];
237  for (unsigned int k = 0; k < n_elems[j]; ++k)
238  {
239  const unsigned int i = i_section_start + k + 1;
240  const Real ref_position = section_start_ref_position + section_dx * k;
241  const Point ref_point = Point(ref_position, 0, 0);
242  element_boundary_coordinates[i] =
243  computeRealPointFromReferencePoint(ref_point, position, R, Rx)(d);
244  }
245 
246  section_start_ref_position += lengths[j];
247  i_section_start += n_elems[j];
248  }
249 
250  return element_boundary_coordinates;
251 }
void mooseError(Args &&... args)
static RealTensorValue computeXRotationTransformationTensor(const Real &rotation)
Computes the rotation transformation tensor.
static RealTensorValue computeDirectionTransformationTensor(const RealVectorValue &dir)
Computes the direction transformation tensor.
static const std::string axis
Definition: NS.h:27
static const std::string R
Definition: NS.h:162
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
MooseEnum getAlignmentAxis() const
Gets an axis MooseEnum for the axis the component is aligned with.
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
Point computeRealPointFromReferencePoint(const Point &p) const
Computes point in 3-D space from a point in reference space.
static const std::string k
Definition: NS.h:130

◆ getLength()

virtual Real DiscreteLineSegmentInterface::getLength ( ) const
inlinevirtual

◆ getNumElems()

virtual Real DiscreteLineSegmentInterface::getNumElems ( ) const
inlinevirtual

Definition at line 31 of file DiscreteLineSegmentInterface.h.

Referenced by HeatTransferFromHeatStructure1Phase::check(), and HeatTransferFromHeatStructure3D1Phase::init().

31 { return _n_elem; }
const unsigned int _n_elem
Total number of axial elements.

◆ getPosition()

virtual Point DiscreteLineSegmentInterface::getPosition ( ) const
inlinevirtual

◆ getRotation()

virtual Real DiscreteLineSegmentInterface::getRotation ( ) const
inlinevirtual

Definition at line 29 of file DiscreteLineSegmentInterface.h.

29 { return _rotation; }
const Real & _rotation
Angle of rotation about the x-axis.

◆ initializeDirectionVector()

RealVectorValue DiscreteLineSegmentInterface::initializeDirectionVector ( const RealVectorValue dir_unnormalized)
staticprivate

Computes a normalized direction vector or reports an error if the zero vector is provided.

Parameters
[in]dir_unnormalizedUnnormalized direction vector

Definition at line 73 of file DiscreteLineSegmentInterface.C.

74 {
75  if (!MooseUtils::absoluteFuzzyEqual(dir_unnormalized.norm(), 0.0))
76  return dir_unnormalized / dir_unnormalized.norm();
77  else
78  mooseError("The parameter 'orientation' must not be the zero vector.");
79 }
auto norm() const -> decltype(std::norm(Real()))
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
void mooseError(Args &&... args)

◆ validParams()

InputParameters DiscreteLineSegmentInterface::validParams ( )
static

Definition at line 17 of file DiscreteLineSegmentInterface.C.

Referenced by GeneratedMeshComponent::validParams(), and DiscreteLineSegmentInterfaceTestAux::validParams().

18 {
20 
21  params.addRequiredParam<Point>("position", "Start position of axis in 3-D space [m]");
23  "orientation",
24  "Direction of axis from start position to end position (no need to normalize)");
25  params.addParam<Real>("rotation", 0.0, "Angle of rotation about the x-axis [degrees]");
26  params.addRequiredParam<std::vector<Real>>("length", "Length of each axial section [m]");
27  params.addRequiredParam<std::vector<unsigned int>>("n_elems",
28  "Number of elements in each axial section");
29 
30  return params;
31 }
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)
InputParameters emptyInputParameters()
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

Member Data Documentation

◆ _dir

const RealVectorValue DiscreteLineSegmentInterface::_dir
protected

Normalized direction of axis from start position to end position.

Definition at line 95 of file DiscreteLineSegmentInterface.h.

Referenced by Component1D::buildMesh(), computeAxialCoordinate(), computeRadialCoordinate(), getAlignmentAxis(), getDirection(), and getElementBoundaryCoordinates().

◆ _dir_unnormalized

const RealVectorValue& DiscreteLineSegmentInterface::_dir_unnormalized
protected

Unnormalized direction of axis from start position to end position.

Definition at line 93 of file DiscreteLineSegmentInterface.h.

◆ _length

Real DiscreteLineSegmentInterface::_length
protected

◆ _lengths

std::vector<Real> DiscreteLineSegmentInterface::_lengths
protected

◆ _moose_object_name_dlsi

const std::string DiscreteLineSegmentInterface::_moose_object_name_dlsi
protected

Name of the MOOSE object.

Definition at line 128 of file DiscreteLineSegmentInterface.h.

Referenced by computeAxialCoordinate().

◆ _n_elem

const unsigned int DiscreteLineSegmentInterface::_n_elem
protected

◆ _n_elems

const std::vector<unsigned int>& DiscreteLineSegmentInterface::_n_elems
protected

◆ _n_sections

const unsigned int DiscreteLineSegmentInterface::_n_sections
protected

◆ _position

const Point& DiscreteLineSegmentInterface::_position
protected

◆ _R

const RealTensorValue DiscreteLineSegmentInterface::_R
protected

Direction transformation tensor.

Definition at line 118 of file DiscreteLineSegmentInterface.h.

Referenced by computeRealPointFromReferencePoint().

◆ _R_inv

const RealTensorValue DiscreteLineSegmentInterface::_R_inv
protected

Inverse direction transformation tensor.

Definition at line 123 of file DiscreteLineSegmentInterface.h.

Referenced by computeReferencePointFromRealPoint().

◆ _rotation

const Real& DiscreteLineSegmentInterface::_rotation
protected

Angle of rotation about the x-axis.

Definition at line 97 of file DiscreteLineSegmentInterface.h.

Referenced by getElementBoundaryCoordinates(), and getRotation().

◆ _Rx

const RealTensorValue DiscreteLineSegmentInterface::_Rx
protected

Rotational transformation tensor about x-axis.

Definition at line 120 of file DiscreteLineSegmentInterface.h.

Referenced by computeRealPointFromReferencePoint().

◆ _Rx_inv

const RealTensorValue DiscreteLineSegmentInterface::_Rx_inv
protected

Inverse rotational transformation tensor about x-axis.

Definition at line 125 of file DiscreteLineSegmentInterface.h.

Referenced by computeReferencePointFromRealPoint().

◆ _section_end

std::vector<Real> DiscreteLineSegmentInterface::_section_end
protected

Axial coordinate of the end of each axial section using the line 'position' as the origin.

Definition at line 112 of file DiscreteLineSegmentInterface.h.

Referenced by DiscreteLineSegmentInterface(), and getAxialSectionIndex().

◆ _x_centers

std::vector<Real> DiscreteLineSegmentInterface::_x_centers
protected

Center axial coordinate of each axial element.

Definition at line 115 of file DiscreteLineSegmentInterface.h.

Referenced by DiscreteLineSegmentInterface(), and getAxialElementIndex().


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