KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
rotate_region_process.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 //
8 // License: BSD License
9 // Kratos default license: kratos/license.txt
10 //
11 // Authors: Aditya Ghantasala, https://github.com/adityaghantasala
12 // Navaneeth K Narayanan
13 //
14 
15 #ifndef ROTATE_REGION_PROCESS_H
16 #define ROTATE_REGION_PROCESS_H
17 
18 // System includes
19 #include <cmath>
20 #include <iostream>
21 #include <string>
22 
23 // External includes
24 
25 // Project includes
26 #include "includes/define.h"
27 #include "containers/model.h"
29 #include "includes/model_part.h"
30 #include "includes/variables.h"
31 #include "processes/process.h"
32 
33 // Application includes
35 
36 namespace Kratos {
37 
48 class KRATOS_API(CHIMERA_APPLICATION) RotateRegionProcess : public Process {
49 public:
52 
55  typedef ProcessInfo::Pointer ProcessInfoPointerType;
56  typedef Matrix MatrixType;
57  typedef Vector VectorType;
58 
60  RotateRegionProcess(ModelPart &rModelPart, Parameters rParameters);
61 
63  ~RotateRegionProcess() = default;
64 
65  void SetAngularVelocity(const double NewAngularVelocity);
66 
67  void ExecuteInitializeSolutionStep() override;
68 
69  std::string Info() const override;
70 
72  void PrintInfo(std::ostream &rOStream) const override;
73 
75  void PrintData(std::ostream& rOStream) const override;
76 
77 protected:
80 
84 
86 
87 private:
90 
94 
104  class RotationSystem {
105  public:
107  KRATOS_CLASS_POINTER_DEFINITION(RotationSystem);
108  /*
109  * @brief Constructor
110  * @param MomentOfInertia The moment of inertia of the system.
111  * @param DampingCoefficient The Damping Coefficient of the system.
112  */
113  RotationSystem(const double MomentOfInertia,
114  const double DampingCoefficient = 0.0);
115 
116  /*
117  * @brief Advances the rotation system in time
118  * @param NewTime Time where the system is to be set
119  */
120  void CloneTimeStep(const double NewTime, const double Dt);
121 
122  /*
123  * @brief Applies the given torque on the system
124  * @param Torque The torque to be applied
125  */
126  void inline ApplyTorque(const double Torque);
127 
128  /*
129  * @brief Calculates the current state of the system
130  */
131  double CalculateCurrentRotationState();
132  /*
133  * @brief Gives the current angle of rotation Theta
134  * @out Current angle of rotation Theta
135  */
136  double GetCurrentTheta() const;
137 
138  /*
139  * @brief Gives the current angular velocity Omega
140  * @out Current angular velocity Omega
141  */
142  double GetCurrentOmega() const;
143 
144  private:
145  double mDt;
146  double mMomentOfInertia;
147  double mDampingCoeff;
148  double mTorque;
149  double mTime;
150  DenseVector<double> mBdf2Coeff;
151  DenseVector<double> mTheta;
152  DenseVector<double> mOmega;
153 
154  /*
155  * @brief Calculates the Intertial torque acting on the system
156  * @out The inertial torque
157  */
158  double CalculateInertiaTorque() const;
159 
160  /*
161  * @brief Calculates the Damping torque acting on the system
162  * @out The damping torque
163  */
164  double CalculateDampingTorque() const;
165 
166  /*
167  * @brief Computes the LHS of the Euler's equations
168  * @out The LHS
169  */
170  double ComputeLHS() const;
171  /*
172  * @brief Computes the RHS of the Euler's equations
173  * @out The RHS
174  */
175  double ComputeRHS() const;
176 
177  /*
178  * @brief Predicts the next time steps theta
179  * @out The predicted value
180  */
181  void Predict();
182 
183  /*
184  * @brief Updates the state of the system with a given update
185  * @param UpdateTheta The update of theta
186  */
187  void Update(double UpdateTheta);
188  };
189 
190  ModelPart &mrModelPart;
191  Parameters mParameters;
192  double mAngularVelocityRadians;
193  array_1d<double, 3> mAxisOfRotationVector;
194  array_1d<double, 3> mCenterOfRotation;
195  double mTheta;
196  double mDTheta;
197  bool mToCalculateTorque;
198  RotationSystem::Pointer mpRotationSystem;
199  double mTime;
200 
202 
205  {
206  return *this;
207  }
208 
209  /*
210  * @brief Calculates the current rotation of modelpart.
211  */
212  void CalculateCurrentRotationState();
213 
214  /*
215  * @brief Calculates the linear velocity v = r x w
216  * @param rAxisOfRotationVector the axis of rotation vector
217  * @param rRadius the radius vector of the point for which v is to be
218  * calculated.
219  * @out rLinearVelocity the calculated linear velocity.
220  */
221  void CalculateLinearVelocity(const DenseVector<double> &rAxisOfRotationVector,
222  const DenseVector<double> &rRadius,
223  DenseVector<double> &rLinearVelocity);
224 
225  /*
226  * @brief Rotates the given node by mTheta around the mAxisOfRotationVector
227  * and mCenterOfRotation This function uses Quaternion for rotations.
228  * @param rCoordinates The nodal coordinates which are to be rotated.
229  * @out rTransformedCoordinates the rotated nodal coordinates.
230  */
231  void TransformNode(const array_1d<double, 3> &rCoordinates,
232  array_1d<double, 3> &rTransformedCoordinates,
233  double Theta) const;
234  /*
235  * @brief Calculates Torque on the given modelpart
236  * @out Torque on the modelpart about the axis of rotation.
237  */
238  double CalculateTorque() const;
239 }; // Class MoveRotorProcess
240 
241 }; // namespace Kratos.
242 
243 #endif // KRATOS_MOVE_ROTOR_PROCESS_H
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
The base class for all processes in Kratos.
Definition: process.h:49
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
This is the to apply rotation to a given modelpart.
Definition: rotate_region_process.h:48
ModelPart::NodeIterator NodeIteratorType
Definition: rotate_region_process.h:53
KRATOS_CLASS_POINTER_DEFINITION(RotateRegionProcess)
Pointer definition of MoveRotorProcess.
Matrix MatrixType
Definition: rotate_region_process.h:56
ProcessInfo ProcessInfoType
Definition: rotate_region_process.h:54
ProcessInfo::Pointer ProcessInfoPointerType
Definition: rotate_region_process.h:55
~RotateRegionProcess()=default
Destructor.
Vector VectorType
Definition: rotate_region_process.h:57
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Dt
Definition: face_heat.py:78
#define KRATOS_CLASS_POINTER_DEFINITION(a)
Definition: smart_pointers.h:69