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

#include <MooseAppCoordTransform.h>

Public Types

enum  Direction : unsigned int { X = 0, Y, Z, INVALID }
 A class scope enumeration for conveniently denoting X, Y, and Z axis directions. More...
 
typedef std::tuple< short int, Real, short int, std::array< Real, 3 >, int, unsigned int, unsigned int, short int, short int, short intMinimalData
 A typedef for conveniency that describes the minimal data necessary to broadcast and build a MooseAppCoordTransform. More...
 

Public Member Functions

 MooseAppCoordTransform ()
 Default constructor. More...
 
 MooseAppCoordTransform (const MooseAppCoordTransform &other)
 
 MooseAppCoordTransform (MooseAppCoordTransform &&other)=default
 
 MooseAppCoordTransform (const MinimalData &minimal_data)
 Construct a coordinate transformation object from the minimal set of data required. More...
 
 MooseAppCoordTransform (const MooseMesh &mesh)
 Construct this object from the provided mesh and its input parameters. More...
 
 ~MooseAppCoordTransform ()=default
 
MooseAppCoordTransformoperator= (const MooseAppCoordTransform &other)
 
MooseAppCoordTransformoperator= (MooseAppCoordTransform &&other)=default
 
MinimalData minimalDataDescription () const
 
Moose::CoordinateSystemType coordinateSystem () const
 
void setTranslationVector (const libMesh::Point &translation)
 Set how much our domain should be translated in order to match a reference frame. More...
 
void setUpDirection (Direction up_direction)
 Will setup a rotation transformation. More...
 
void setRotation (Real alpha, Real beta, Real gamma)
 Setup an extrinsic rotation defined in the following way: More...
 
void setLengthUnit (const MooseUnits &length_unit)
 Set the scaling transformation. More...
 
const MooseUnitslengthUnit () const
 
void setCoordinateSystem (Moose::CoordinateSystemType system_type, Direction rz_symmetry_axis=INVALID)
 Set our coordinate system. More...
 
void setCoordinateSystem (const MooseMesh &mesh)
 Set our coordinate system based on the MooseMesh coordinate system data. More...
 
void computeRS ()
 Compute the RS and (RS)^{-1} matrices. More...
 
void transformMesh (MooseMesh &mesh, const libMesh::Point &translation)
 Transforms the entire mesh with the coordinate transform This can be done to output in position, or to avoid transforming on every data point. More...
 
bool hasScalingOrRotationTransformation () const
 Returns true if the app has scaling and/or rotation transformation. More...
 

Static Public Member Functions

static InputParameters validParams ()
 Describes the parameters this object can take to setup transformations. More...
 

Private Member Functions

Direction processZAxis (Direction z_axis)
 If the coordinate system type is RZ, then we return the provided argument. More...
 

Private Attributes

std::unique_ptr< libMesh::RealTensorValue_scale
 Represents a forward scaling transformation from our units to reference frame units of meters. More...
 
std::unique_ptr< libMesh::RealTensorValue_rotate
 Represents a forward rotation transformation from our domain to the reference frame domain. More...
 
std::unique_ptr< libMesh::RealTensorValue_rs
 Represents the product of rotation and scaling transformations. More...
 
std::unique_ptr< libMesh::RealTensorValue_rs_inverse
 Represents the inverse of the product of rotation and scaling transformations. More...
 
Moose::CoordinateSystemType _coord_type
 Our coordinate system. More...
 
Direction _r_axis
 If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate. More...
 
Direction _z_axis
 If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate. More...
 
bool _has_different_coord_sys
 Whether we have different coordinate systems within our single domain. More...
 
bool _using_general_rz_coord_axes
 Whether general axisymmetric coordinate axes are being used. More...
 
MooseUnits _length_unit
 How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches. More...
 
std::array< Real, 3 > _euler_angles
 The Euler angles describing rotation. More...
 
bool _mesh_transformed
 Whether the mesh has been translated and rotated. More...
 

Friends

class MultiAppCoordTransform
 

Detailed Description

Definition at line 24 of file MooseAppCoordTransform.h.

Member Typedef Documentation

◆ MinimalData

typedef std::tuple<short int, Real, short int, std::array<Real, 3>, int, unsigned int, unsigned int, short int, short int, short int> MooseAppCoordTransform::MinimalData

A typedef for conveniency that describes the minimal data necessary to broadcast and build a MooseAppCoordTransform.

The data is such: 0: whether a scaling matrix exists 1: a value describing the scaling 2: whether a rotation matrix exists 3: the Euler angles describing the rotation 4: the coordinate system type 5: the r-axis direction 6: the z-axis direction 7: whether there are multiple coordinate system types on the mesh 8: whether general axisymmetric coordinate axes are being used 9: whether the mesh has been transformed using the transform

Definition at line 62 of file MooseAppCoordTransform.h.

Member Enumeration Documentation

◆ Direction

A class scope enumeration for conveniently denoting X, Y, and Z axis directions.

Enumerator
INVALID 

Definition at line 30 of file MooseAppCoordTransform.h.

Constructor & Destructor Documentation

◆ MooseAppCoordTransform() [1/5]

MooseAppCoordTransform::MooseAppCoordTransform ( )

Default constructor.

If no other methods are called to set rotation, translation, or scaling, then when operator() is called the result will be the passed-in point, e.g. no transformation will occur

Definition at line 305 of file MooseAppCoordTransform.C.

307  _r_axis(INVALID),
308  _z_axis(INVALID),
311  _length_unit(std::string("1*m")),
312  _euler_angles(),
313  _mesh_transformed(false)
314 {
315 }
bool _mesh_transformed
Whether the mesh has been translated and rotated.
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
Moose::CoordinateSystemType _coord_type
Our coordinate system.
bool _has_different_coord_sys
Whether we have different coordinate systems within our single domain.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
MooseUnits _length_unit
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ MooseAppCoordTransform() [2/5]

MooseAppCoordTransform::MooseAppCoordTransform ( const MooseAppCoordTransform other)

Definition at line 317 of file MooseAppCoordTransform.C.

318  : _coord_type(other._coord_type),
319  _r_axis(other._r_axis),
320  _z_axis(other._z_axis),
323  _length_unit(other._length_unit),
326 {
327  if (other._scale)
328  _scale = std::make_unique<RealTensorValue>(*other._scale);
329  if (other._rotate)
330  _rotate = std::make_unique<RealTensorValue>(*other._rotate);
331  computeRS();
332 }
std::unique_ptr< libMesh::RealTensorValue > _scale
Represents a forward scaling transformation from our units to reference frame units of meters...
bool _mesh_transformed
Whether the mesh has been translated and rotated.
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
Moose::CoordinateSystemType _coord_type
Our coordinate system.
bool _has_different_coord_sys
Whether we have different coordinate systems within our single domain.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...
MooseUnits _length_unit
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.
void computeRS()
Compute the RS and (RS)^{-1} matrices.
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ MooseAppCoordTransform() [3/5]

MooseAppCoordTransform::MooseAppCoordTransform ( MooseAppCoordTransform &&  other)
default

◆ MooseAppCoordTransform() [4/5]

MooseAppCoordTransform::MooseAppCoordTransform ( const MinimalData minimal_data)

Construct a coordinate transformation object from the minimal set of data required.

Definition at line 334 of file MooseAppCoordTransform.C.

335  : _coord_type(static_cast<Moose::CoordinateSystemType>(std::get<4>(minimal_data))),
336  _r_axis(static_cast<Direction>(std::get<5>(minimal_data))),
337  _z_axis(static_cast<Direction>(std::get<6>(minimal_data))),
338  _has_different_coord_sys(std::get<7>(minimal_data)),
339  _using_general_rz_coord_axes(std::get<8>(minimal_data)),
340  _length_unit(std::string("1*m")),
341  _euler_angles(std::get<3>(minimal_data)),
342  _mesh_transformed(std::get<9>(minimal_data))
343 {
344  if (std::get<0>(minimal_data))
345  setLengthUnit(MooseUnits(std::to_string(std::get<1>(minimal_data)) + "*m"));
346  if (std::get<2>(minimal_data))
347  _rotate = std::make_unique<RealTensorValue>(RealTensorValue::extrinsic_rotation_matrix(
349  computeRS();
350 }
static TensorValue< Real > extrinsic_rotation_matrix(Real angle1_deg, Real angle2_deg, Real angle3_deg)
bool _mesh_transformed
Whether the mesh has been translated and rotated.
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
Moose::CoordinateSystemType _coord_type
Our coordinate system.
bool _has_different_coord_sys
Whether we have different coordinate systems within our single domain.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...
void setLengthUnit(const MooseUnits &length_unit)
Set the scaling transformation.
MooseUnits _length_unit
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.
Physical unit management class with runtime unit string parsing, unit checking, unit conversion...
Definition: Units.h:32
void computeRS()
Compute the RS and (RS)^{-1} matrices.
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ MooseAppCoordTransform() [5/5]

MooseAppCoordTransform::MooseAppCoordTransform ( const MooseMesh mesh)

Construct this object from the provided mesh and its input parameters.

See the validParams implementation for valid parameters

Definition at line 259 of file MooseAppCoordTransform.C.

261  _r_axis(INVALID),
262  _z_axis(INVALID),
265  _length_unit(std::string("1*m")),
266  _euler_angles(),
267  _mesh_transformed(false)
268 {
269  //
270  // Coordinate system transformation
271  //
272  setCoordinateSystem(mesh);
273 
274  const auto & params = mesh.parameters();
275 
276  //
277  // Rotation
278  //
279  const bool has_alpha = params.isParamValid("alpha_rotation");
280  const bool has_beta = params.isParamValid("beta_rotation");
281  const bool has_gamma = params.isParamValid("gamma_rotation");
282  const auto & up_direction = params.get<MooseEnum>("up_direction");
283 
284  if (has_alpha || has_beta || has_gamma)
285  {
286  if (up_direction.isValid())
287  mooseError("Cannot simultaneously set rotation angles as well as an up direction");
288 
289  const auto alpha = (has_alpha ? params.get<Real>("alpha_rotation") : Real(0));
290  const auto beta = (has_beta ? params.get<Real>("beta_rotation") : Real(0));
291  const auto gamma = (has_gamma ? params.get<Real>("gamma_rotation") : Real(0));
292 
293  setRotation(alpha, beta, gamma);
294  }
295  else if (up_direction.isValid())
296  setUpDirection(Direction(static_cast<unsigned int>(int(up_direction))));
297 
298  //
299  // Scaling
300  //
301  if (params.isParamValid("length_unit"))
302  setLengthUnit(MooseUnits(params.get<std::string>("length_unit")));
303 }
bool _mesh_transformed
Whether the mesh has been translated and rotated.
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
Moose::CoordinateSystemType _coord_type
Our coordinate system.
MeshBase & mesh
Direction
A class scope enumeration for conveniently denoting X, Y, and Z axis directions.
bool _has_different_coord_sys
Whether we have different coordinate systems within our single domain.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
void setLengthUnit(const MooseUnits &length_unit)
Set the scaling transformation.
void setCoordinateSystem(Moose::CoordinateSystemType system_type, Direction rz_symmetry_axis=INVALID)
Set our coordinate system.
MooseUnits _length_unit
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.
Physical unit management class with runtime unit string parsing, unit checking, unit conversion...
Definition: Units.h:32
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void setRotation(Real alpha, Real beta, Real gamma)
Setup an extrinsic rotation defined in the following way:
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.
void setUpDirection(Direction up_direction)
Will setup a rotation transformation.

◆ ~MooseAppCoordTransform()

MooseAppCoordTransform::~MooseAppCoordTransform ( )
default

Member Function Documentation

◆ computeRS()

void MooseAppCoordTransform::computeRS ( )

Compute the RS and (RS)^{-1} matrices.

Definition at line 395 of file MooseAppCoordTransform.C.

Referenced by MooseAppCoordTransform(), and operator=().

396 {
397  if (_scale || _rotate)
398  {
399  _rs = std::make_unique<RealTensorValue>(RealTensorValue(1, 0, 0, 0, 1, 0, 0, 0, 1));
400 
401  if (_scale)
402  *_rs = *_scale * *_rs;
403  if (_rotate)
404  *_rs = *_rotate * *_rs;
405 
406  _rs_inverse = std::make_unique<RealTensorValue>(_rs->inverse());
407  }
408  else
409  {
410  _rs.reset();
411  _rs_inverse.reset();
412  }
413 }
std::unique_ptr< libMesh::RealTensorValue > _scale
Represents a forward scaling transformation from our units to reference frame units of meters...
std::unique_ptr< libMesh::RealTensorValue > _rs_inverse
Represents the inverse of the product of rotation and scaling transformations.
std::unique_ptr< libMesh::RealTensorValue > _rs
Represents the product of rotation and scaling transformations.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...

◆ coordinateSystem()

Moose::CoordinateSystemType MooseAppCoordTransform::coordinateSystem ( ) const
inline
Returns
our coordinate system

Definition at line 107 of file MooseAppCoordTransform.h.

107 { return _coord_type; }
Moose::CoordinateSystemType _coord_type
Our coordinate system.

◆ hasScalingOrRotationTransformation()

bool MooseAppCoordTransform::hasScalingOrRotationTransformation ( ) const

Returns true if the app has scaling and/or rotation transformation.

Definition at line 443 of file MooseAppCoordTransform.C.

Referenced by FEProblemBase::checkProblemIntegrity(), MultiAppCoordTransform::hasNonTranslationTransformation(), MultiAppCoordTransform::setDestinationCoordTransform(), and transformMesh().

444 {
445  if (_rs)
446  for (const auto i : make_range(Moose::dim))
447  for (const auto j : make_range(Moose::dim))
448  {
449  const auto matrix_elem = (*_rs)(i, j);
450  if (i == j)
451  {
452  if (!MooseUtils::absoluteFuzzyEqual(matrix_elem, 1))
453  return true;
454  }
455  else if (!MooseUtils::absoluteFuzzyEqual(matrix_elem, 0))
456  return true;
457  }
458 
459  return false;
460 }
bool absoluteFuzzyEqual(const T &var1, const T2 &var2, const T3 &tol=libMesh::TOLERANCE *libMesh::TOLERANCE)
Function to check whether two variables are equal within an absolute tolerance.
Definition: MooseUtils.h:380
static constexpr std::size_t dim
This is the dimension of all vector and tensor datastructures used in MOOSE.
Definition: Moose.h:153
std::unique_ptr< libMesh::RealTensorValue > _rs
Represents the product of rotation and scaling transformations.
IntRange< T > make_range(T beg, T end)

◆ lengthUnit()

const MooseUnits& MooseAppCoordTransform::lengthUnit ( ) const
inline
Returns
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches

Definition at line 158 of file MooseAppCoordTransform.h.

158 { return _length_unit; }
MooseUnits _length_unit
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.

◆ minimalDataDescription()

MooseAppCoordTransform::MinimalData MooseAppCoordTransform::minimalDataDescription ( ) const
Returns
the minimal data necessary to describe this coordinate transformation object. This data can be broadcast and used to construct identical coordinate transformation objects on other processes

Definition at line 353 of file MooseAppCoordTransform.C.

354 {
355  const Real scale_factor = _scale ? (*_scale)(0, 0) : 1;
356  return {static_cast<short int>(bool(_scale)),
357  scale_factor,
358  static_cast<short int>(bool(_rotate)),
360  static_cast<int>(_coord_type),
361  static_cast<unsigned int>(_r_axis),
362  static_cast<unsigned int>(_z_axis),
363  static_cast<short int>(_has_different_coord_sys),
364  static_cast<short int>(_using_general_rz_coord_axes),
365  static_cast<short int>(_mesh_transformed)};
366 }
std::unique_ptr< libMesh::RealTensorValue > _scale
Represents a forward scaling transformation from our units to reference frame units of meters...
bool _mesh_transformed
Whether the mesh has been translated and rotated.
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
Moose::CoordinateSystemType _coord_type
Our coordinate system.
bool _has_different_coord_sys
Whether we have different coordinate systems within our single domain.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ operator=() [1/2]

MooseAppCoordTransform & MooseAppCoordTransform::operator= ( const MooseAppCoordTransform other)

Definition at line 369 of file MooseAppCoordTransform.C.

370 {
371  _coord_type = other._coord_type;
372  _r_axis = other._r_axis;
373  _z_axis = other._z_axis;
376  _length_unit = other._length_unit;
379 
380  if (other._scale)
381  _scale = std::make_unique<RealTensorValue>(*other._scale);
382  else
383  _scale.reset();
384  if (other._rotate)
385  _rotate = std::make_unique<RealTensorValue>(*other._rotate);
386  else
387  _rotate.reset();
388 
389  computeRS();
390 
391  return *this;
392 }
std::unique_ptr< libMesh::RealTensorValue > _scale
Represents a forward scaling transformation from our units to reference frame units of meters...
bool _mesh_transformed
Whether the mesh has been translated and rotated.
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
Moose::CoordinateSystemType _coord_type
Our coordinate system.
bool _has_different_coord_sys
Whether we have different coordinate systems within our single domain.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...
MooseUnits _length_unit
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.
void computeRS()
Compute the RS and (RS)^{-1} matrices.
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ operator=() [2/2]

MooseAppCoordTransform& MooseAppCoordTransform::operator= ( MooseAppCoordTransform &&  other)
default

◆ processZAxis()

MooseAppCoordTransform::Direction MooseAppCoordTransform::processZAxis ( Direction  z_axis)
private

If the coordinate system type is RZ, then we return the provided argument.

Otherwise we return INVALID

Definition at line 20 of file MooseAppCoordTransform.C.

21 {
22  return _coord_type == Moose::COORD_RZ ? z_axis : INVALID;
23 }
Moose::CoordinateSystemType _coord_type
Our coordinate system.

◆ setCoordinateSystem() [1/2]

void MooseAppCoordTransform::setCoordinateSystem ( Moose::CoordinateSystemType  system_type,
Direction  rz_symmetry_axis = INVALID 
)

Set our coordinate system.

Parameters
system_typethe coordinate system type
rz_symmetry_axisthe axial coordinate, e.g. the axis of symmetry

Definition at line 143 of file MooseAppCoordTransform.C.

Referenced by MooseAppCoordTransform().

145 {
146  _coord_type = coord_type;
147 
149  {
150  if (rz_symmetry_axis == INVALID)
151  mooseError("For RZ coordinate systems, the 'rz_symmetry_axis' parameter must be provided to "
152  "'MooseAppCoordTransform::setCoordinateSystem'");
153 
154  _z_axis = rz_symmetry_axis;
155  _r_axis = _z_axis == X ? Y : X;
156  }
158  _r_axis = X;
159 }
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
Moose::CoordinateSystemType _coord_type
Our coordinate system.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.

◆ setCoordinateSystem() [2/2]

void MooseAppCoordTransform::setCoordinateSystem ( const MooseMesh mesh)

Set our coordinate system based on the MooseMesh coordinate system data.

Definition at line 162 of file MooseAppCoordTransform.C.

163 {
164  const auto & params = mesh.parameters();
165 
166  // If we have multiple different coordinate system types in our problem, we
167  // take note of it because that can cause issues if there is a non-Cartesian
168  // destination coordinate system
169  const auto & coord_sys = mesh.getCoordSystem();
170  std::unordered_set<Moose::CoordinateSystemType> coord_types;
171  auto map_it = coord_sys.begin();
172  // It's possible that the mesh is not in a complete state
173  if (map_it == coord_sys.end())
175  else
177  map_it->second,
178  Direction(static_cast<unsigned int>(int(params.get<MooseEnum>("rz_coord_axis")))));
179  for (; map_it != coord_sys.end(); ++map_it)
180  coord_types.insert(map_it->second);
181 
182  _has_different_coord_sys = coord_types.size() > 1;
183 
184  if (mesh.usingGeneralAxisymmetricCoordAxes())
186 }
MeshBase & mesh
Direction
A class scope enumeration for conveniently denoting X, Y, and Z axis directions.
bool _has_different_coord_sys
Whether we have different coordinate systems within our single domain.
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
void setCoordinateSystem(Moose::CoordinateSystemType system_type, Direction rz_symmetry_axis=INVALID)
Set our coordinate system.
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.

◆ setLengthUnit()

void MooseAppCoordTransform::setLengthUnit ( const MooseUnits length_unit)

Set the scaling transformation.

Parameters
length_unitHow much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches. We will save off the value provided to this in the _length_unit data member as well as set the scaling transform

Definition at line 133 of file MooseAppCoordTransform.C.

Referenced by MooseAppCoordTransform().

134 {
135  _length_unit = length_unit;
136  const auto scale = Real(_length_unit / MooseUnits("m"));
137  _scale =
138  std::make_unique<RealTensorValue>(RealTensorValue(scale, 0, 0, 0, scale, 0, 0, 0, scale));
139  computeRS();
140 }
std::unique_ptr< libMesh::RealTensorValue > _scale
Represents a forward scaling transformation from our units to reference frame units of meters...
void scale(MeshBase &mesh, const Real xs, const Real ys=0., const Real zs=0.)
MooseUnits _length_unit
How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.
Physical unit management class with runtime unit string parsing, unit checking, unit conversion...
Definition: Units.h:32
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void computeRS()
Compute the RS and (RS)^{-1} matrices.

◆ setRotation()

void MooseAppCoordTransform::setRotation ( Real  alpha,
Real  beta,
Real  gamma 
)

Setup an extrinsic rotation defined in the following way:

  1. rotate by alpha degrees about the z-axis
  2. rotate by beta degrees about the x-axis
  3. rotate by gamma degrees about the z-axis Definitions of the resulting matrix are found in the last row of the Proper Euler angles column of https://en.wikipedia.org/wiki/Euler_angles#Rotation_matrix. These rotation angles should describe how points in our domain should be rotated in order to arrive back in the reference frame. For instance, in 2D your mesh may appear 90 degrees rotated (around the z-axis) with respect to the reference frame. In such a case, the angle set you should provide to this function is {-90, 0, 0}, e.g. provide forward transformation angles that will map points from your domain to the reference domain

If our coordinate system is RZ, then only certain values of alpha, beta, and gamma will be accepted such that the radial and axial coordinates are rotated onto Cartesian axes and the resulting radial coordinate is non-negative

Definition at line 91 of file MooseAppCoordTransform.C.

Referenced by MooseAppCoordTransform().

92 {
93  const bool must_rotate_axes =
95  bool axes_rotated = false;
96  if (must_rotate_axes)
97  {
98  const auto angles = std::make_tuple(alpha, beta, gamma);
99  if (angles == std::make_tuple(0, 90, 0))
100  {
101  if (_r_axis == X)
102  {
103  mooseAssert((_coord_type == Moose::COORD_RZ && _z_axis == Y) ||
105  "'_z_axis' is not an expected value");
106  _r_axis = X;
108  }
109  else if (_r_axis == Y)
110  {
111  mooseAssert((_coord_type == Moose::COORD_RZ && _z_axis == X) ||
113  "'_z_axis' is not an expected value");
114  _r_axis = Z;
116  }
117  axes_rotated = true;
118  }
119  }
120 
121  _euler_angles = {{alpha, beta, gamma}};
122  _rotate = std::make_unique<RealTensorValue>(
123  RealTensorValue::extrinsic_rotation_matrix(alpha, beta, gamma));
124  computeRS();
125 
126  if (must_rotate_axes && !axes_rotated)
127  mooseError("Unsupported manual angle prescription in 'MooseAppCoordTransform::setRotation'. "
128  "For non-Cartesian coordinate systems, the only currently supported rotation is "
129  "(alpha, beta, gamma) = (0, 90, 0)");
130 }
static TensorValue< Real > extrinsic_rotation_matrix(Real angle1_deg, Real angle2_deg, Real angle3_deg)
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
Moose::CoordinateSystemType _coord_type
Our coordinate system.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...
void computeRS()
Compute the RS and (RS)^{-1} matrices.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ setTranslationVector()

void MooseAppCoordTransform::setTranslationVector ( const libMesh::Point translation)

Set how much our domain should be translated in order to match a reference frame.

In practice we choose the parent application to be the reference frame with respect to translation, e.g. the parent application origin is the reference frame origin, and we set the translation vectors of child applications to the multiapp positions parameter. Similarly to the setRotation with angles API, this represents a forward transformation from our domain to the reference domain

◆ setUpDirection()

void MooseAppCoordTransform::setUpDirection ( Direction  up_direction)

Will setup a rotation transformation.

The rotation transformation will be a single 90-degree rotation defined such that a point on the axis specified by up_direction is rotated onto the Y-axis, which is our canonical/reference-frame up-direction

Parameters
up_directionWhat direction corresponds to "up" (e.g. the opposite direction of gravity) in our moose mesh

Definition at line 26 of file MooseAppCoordTransform.C.

Referenced by MooseAppCoordTransform().

27 {
28  Real alpha = 0, beta = 0, gamma = 0;
29 
30  const bool must_rotate_axes =
32  // Don't error immediately for unit testing purposes
33  bool negative_radii = false;
34 
35  if (up_direction == X)
36  {
37  alpha = 90, beta = 0, gamma = 0;
38  if (must_rotate_axes)
39  {
40  if (_r_axis == X)
41  {
42  _r_axis = Y;
44  }
45  else if (_r_axis == Y)
46  {
47  negative_radii = true;
48  _r_axis = X;
50  }
51  else
52  mooseError("Bad r-axis value");
53  }
54  }
55  else if (up_direction == Y)
56  alpha = 0, beta = 0, gamma = 0;
57  else if (up_direction == Z)
58  {
59  alpha = 0, beta = -90, gamma = 0;
60  if (must_rotate_axes)
61  {
62  if (_r_axis == X)
63  {
64  _r_axis = X;
66  }
67  else if (_r_axis == Y)
68  {
69  negative_radii = true;
70  _r_axis = Z;
72  }
73  else
74  mooseError("Bad r-axis value");
75  }
76  }
77  else
78  mooseError("Bad up direction value");
79 
80  _euler_angles = {{alpha, beta, gamma}};
81 
82  _rotate = std::make_unique<RealTensorValue>(
83  RealTensorValue::extrinsic_rotation_matrix(alpha, beta, gamma));
84  computeRS();
85 
86  if (negative_radii)
87  mooseError("Rotation yields negative radial values");
88 }
static TensorValue< Real > extrinsic_rotation_matrix(Real angle1_deg, Real angle2_deg, Real angle3_deg)
Direction _z_axis
If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
Moose::CoordinateSystemType _coord_type
Our coordinate system.
Direction _r_axis
If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.
Direction processZAxis(Direction z_axis)
If the coordinate system type is RZ, then we return the provided argument.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void computeRS()
Compute the RS and (RS)^{-1} matrices.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ transformMesh()

void MooseAppCoordTransform::transformMesh ( MooseMesh mesh,
const libMesh::Point translation 
)

Transforms the entire mesh with the coordinate transform This can be done to output in position, or to avoid transforming on every data point.

Parameters
meshthe mesh to modify, usually the child app mesh
translationthe translation to apply to the mesh, often the app position

Definition at line 416 of file MooseAppCoordTransform.C.

417 {
418  // Transforming a RZ or R-spherical mesh doesnt always make sense, disallow it
420  mooseError("Running MultiApps 'in position' is only supported for XYZ coordinate systems");
421 
423  mooseError("Scaling and rotation are currently not supported for general axisymmetric "
424  "coordinate systems.");
425 
426  if (_mesh_transformed)
427  mooseError("App mesh is being transformed twice");
428 
429  // Apply all the transformation to the mesh
430  if (_scale)
431  MeshTools::Modification::scale(mesh, (*_scale)(0, 0), (*_scale)(1, 1), (*_scale)(2, 2));
432  if (_rotate)
433  MeshTools::Modification::rotate(mesh, _euler_angles[0], _euler_angles[1], _euler_angles[2]);
434  if (translation != Point(0, 0, 0))
435  MeshTools::Modification::translate(mesh, translation(0), translation(1), translation(2));
436 
437  // Translation, scaling and rotation need not be applied anymore when performing coordinate
438  // transforms
439  _mesh_transformed = true;
440 }
std::unique_ptr< libMesh::RealTensorValue > _scale
Represents a forward scaling transformation from our units to reference frame units of meters...
bool _mesh_transformed
Whether the mesh has been translated and rotated.
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:302
Moose::CoordinateSystemType _coord_type
Our coordinate system.
bool hasScalingOrRotationTransformation() const
Returns true if the app has scaling and/or rotation transformation.
std::unique_ptr< libMesh::RealTensorValue > _rotate
Represents a forward rotation transformation from our domain to the reference frame domain...
bool _using_general_rz_coord_axes
Whether general axisymmetric coordinate axes are being used.
std::array< Real, 3 > _euler_angles
The Euler angles describing rotation.

◆ validParams()

InputParameters MooseAppCoordTransform::validParams ( )
static

Describes the parameters this object can take to setup transformations.

These include parameters related to coordinate system type, rotation, and scaling

One entry of coord system per block, the size of _blocks and _coord_sys has to match, except:

  1. _blocks.size() == 0, then there needs to be just one entry in _coord_sys, which will be set for the whole domain
  2. _blocks.size() > 0 and no coordinate system was specified, then the whole domain will be XYZ.
  3. _blocks.size() > 0 and one coordinate system was specified, then the whole domain will be that system.

Definition at line 189 of file MooseAppCoordTransform.C.

Referenced by FEProblemBase::checkProblemIntegrity(), and MooseMesh::validParams().

190 {
191  auto params = emptyInputParameters();
197  params.addDeprecatedParam<std::vector<SubdomainName>>(
198  "block",
199  "Block IDs for the coordinate systems.",
200  "Please use the 'coord_block' parameter instead.");
201  params.addParam<std::vector<SubdomainName>>(
202  "coord_block",
203  "Block IDs for the coordinate systems. If this parameter is specified, then it must "
204  "encompass all the subdomains on the mesh.");
205  MultiMooseEnum coord_types("XYZ RZ RSPHERICAL", "XYZ");
206  MooseEnum rz_coord_axis("X=0 Y=1", "Y");
207  params.addParam<MultiMooseEnum>(
208  "coord_type", coord_types, "Type of the coordinate system per block param");
209  params.addParam<MooseEnum>(
210  "rz_coord_axis", rz_coord_axis, "The rotation axis (X | Y) for axisymmetric coordinates");
211  params.addParam<std::vector<SubdomainName>>(
212  "rz_coord_blocks", "Blocks using general axisymmetric coordinate systems");
213  params.addParam<std::vector<Point>>("rz_coord_origins",
214  "Axis origin points for each block in 'rz_coord_blocks'");
215  params.addParam<std::vector<RealVectorValue>>(
216  "rz_coord_directions", "Axis directions for each block in 'rz_coord_blocks'");
217  params.addParam<std::string>(
218  "length_unit",
219  "How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5inches");
220  params.addRangeCheckedParam<Real>(
221  "alpha_rotation",
222  "-180<alpha_rotation<=180",
223  "The number of degrees that the domain should be alpha-rotated using the Euler "
224  "angle ZXZ convention from https://en.wikipedia.org/wiki/Euler_angles#Rotation_matrix in "
225  "order to align with a canonical physical space of your choosing.");
226  params.addRangeCheckedParam<Real>(
227  "beta_rotation",
228  "-180<beta_rotation<=180",
229  "The number of degrees that the domain should be beta-rotated using the Euler "
230  "angle ZXZ convention from https://en.wikipedia.org/wiki/Euler_angles#Rotation_matrix in "
231  "order to align with a canonical physical space of your choosing.");
232  params.addRangeCheckedParam<Real>(
233  "gamma_rotation",
234  "-180<gamma_rotation<=180",
235  "The number of degrees that the domain should be gamma-rotated using the Euler "
236  "angle ZXZ convention from https://en.wikipedia.org/wiki/Euler_angles#Rotation_matrix in "
237  "order to align with a canonical physical space of your choosing.");
238  MooseEnum up_direction("X=0 Y=1 Z=2");
239  params.addParam<MooseEnum>(
240  "up_direction",
241  up_direction,
242  "Specify what axis corresponds to the up direction in physical space (the opposite of the "
243  "gravity vector if you will). If this parameter is provided, we will perform a single 90 "
244  "degree rotation of the domain--if the provided axis is 'x' or 'z', we will not rotate if "
245  "the axis is 'y'--such that a point which was on the provided axis will now lie on the "
246  "y-axis, e.g. the y-axis is our canonical up direction. If you want finer grained control "
247  "than this, please use the 'alpha_rotation', 'beta_rotation', and 'gamma_rotation' "
248  "parameters.");
249  params.addParamNamesToGroup(
250  "block coord_type coord_block rz_coord_axis rz_coord_blocks rz_coord_origins "
251  "rz_coord_directions",
252  "Coordinate system");
253  params.addParamNamesToGroup(
254  "length_unit alpha_rotation beta_rotation gamma_rotation up_direction",
255  "Transformations relative to parent application frame of reference");
256  return params;
257 }
InputParameters emptyInputParameters()
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:33
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type...

Friends And Related Function Documentation

◆ MultiAppCoordTransform

friend class MultiAppCoordTransform
friend

Definition at line 238 of file MooseAppCoordTransform.h.

Member Data Documentation

◆ _coord_type

Moose::CoordinateSystemType MooseAppCoordTransform::_coord_type
private

◆ _euler_angles

std::array<Real, 3> MooseAppCoordTransform::_euler_angles
private

The Euler angles describing rotation.

Definition at line 232 of file MooseAppCoordTransform.h.

Referenced by minimalDataDescription(), MooseAppCoordTransform(), operator=(), and transformMesh().

◆ _has_different_coord_sys

bool MooseAppCoordTransform::_has_different_coord_sys
private

Whether we have different coordinate systems within our single domain.

If we do, this will be problematic if we need to collapse from our space into an RZ or RSPHERICAL space because we are only ever provided with a point argument and not a subdomain ID argument. Consequently we will not know in what coordinate system our point lies and will not know how to perform the dimension collapse, and so we will error

Definition at line 223 of file MooseAppCoordTransform.h.

Referenced by minimalDataDescription(), operator=(), and MultiAppCoordTransform::setDestinationCoordTransform().

◆ _length_unit

MooseUnits MooseAppCoordTransform::_length_unit
private

How much distance one mesh length unit represents, e.g. 1 cm, 1 nm, 1 ft, 5 inches.

Definition at line 229 of file MooseAppCoordTransform.h.

Referenced by lengthUnit(), and operator=().

◆ _mesh_transformed

bool MooseAppCoordTransform::_mesh_transformed
private

Whether the mesh has been translated and rotated.

In this case, applying the transform every time is no longer necessary

Definition at line 236 of file MooseAppCoordTransform.h.

Referenced by MultiAppCoordTransform::mapBack(), minimalDataDescription(), MultiAppCoordTransform::operator()(), operator=(), and transformMesh().

◆ _r_axis

Direction MooseAppCoordTransform::_r_axis
private

If we are RZ or RSPHERICAL, the Cartesian axis corresponding to the radial coordinate.

Definition at line 214 of file MooseAppCoordTransform.h.

Referenced by minimalDataDescription(), MultiAppCoordTransform::operator()(), and operator=().

◆ _rotate

std::unique_ptr<libMesh::RealTensorValue> MooseAppCoordTransform::_rotate
private

Represents a forward rotation transformation from our domain to the reference frame domain.

Definition at line 203 of file MooseAppCoordTransform.h.

Referenced by computeRS(), minimalDataDescription(), MooseAppCoordTransform(), operator=(), and transformMesh().

◆ _rs

std::unique_ptr<libMesh::RealTensorValue> MooseAppCoordTransform::_rs
private

Represents the product of rotation and scaling transformations.

Definition at line 206 of file MooseAppCoordTransform.h.

Referenced by computeRS(), hasScalingOrRotationTransformation(), and MultiAppCoordTransform::operator()().

◆ _rs_inverse

std::unique_ptr<libMesh::RealTensorValue> MooseAppCoordTransform::_rs_inverse
private

Represents the inverse of the product of rotation and scaling transformations.

Definition at line 209 of file MooseAppCoordTransform.h.

Referenced by computeRS(), and MultiAppCoordTransform::mapBack().

◆ _scale

std::unique_ptr<libMesh::RealTensorValue> MooseAppCoordTransform::_scale
private

Represents a forward scaling transformation from our units to reference frame units of meters.

This matrix will be diagonal

Definition at line 200 of file MooseAppCoordTransform.h.

Referenced by computeRS(), minimalDataDescription(), MooseAppCoordTransform(), operator=(), and transformMesh().

◆ _using_general_rz_coord_axes

bool MooseAppCoordTransform::_using_general_rz_coord_axes
private

Whether general axisymmetric coordinate axes are being used.

Definition at line 226 of file MooseAppCoordTransform.h.

Referenced by minimalDataDescription(), operator=(), MultiAppCoordTransform::setDestinationCoordTransform(), and transformMesh().

◆ _z_axis

Direction MooseAppCoordTransform::_z_axis
private

If we are RZ, the Cartesian axis corresponding to the axial/axis-of-symmetry coordinate.

Definition at line 216 of file MooseAppCoordTransform.h.

Referenced by minimalDataDescription(), MultiAppCoordTransform::operator()(), and operator=().


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