ross.Rotor
Contents
ross.Rotor#
- class ross.Rotor(shaft_elements, disk_elements=None, bearing_elements=None, point_mass_elements=None, min_w=None, max_w=None, rated_w=None, tag=None)#
A rotor object.
This class will create a rotor with the shaft, disk, bearing and seal elements provided.
- Parameters:
- shaft_elementslist
List with the shaft elements
- disk_elementslist
List with the disk elements
- bearing_elementslist
List with the bearing elements
- point_mass_elements: list
List with the point mass elements
- tagstr
A tag for the rotor
- Attributes:
- evaluesarray
Rotor’s eigenvalues.
- evectorsarray
Rotor’s eigenvectors.
- wnarray
Rotor’s natural frequencies in rad/s.
- wdarray
Rotor’s damped natural frequencies in rad/s.
- Returns:
- A rotor object.
Examples
>>> # Rotor without damping with 2 shaft elements 1 disk and 2 bearings >>> import ross as rs >>> steel = rs.materials.steel >>> z = 0 >>> le = 0.25 >>> i_d = 0 >>> o_d = 0.05 >>> tim0 = rs.ShaftElement(le, i_d, o_d, ... material=steel, ... shear_effects=True, ... rotary_inertia=True, ... gyroscopic=True) >>> tim1 = rs.ShaftElement(le, i_d, o_d, ... material=steel, ... shear_effects=True, ... rotary_inertia=True, ... gyroscopic=True) >>> shaft_elm = [tim0, tim1] >>> disk0 = rs.DiskElement.from_geometry(1, steel, 0.07, 0.05, 0.28) >>> stf = 1e6 >>> bearing0 = rs.BearingElement(0, kxx=stf, cxx=0) >>> bearing1 = rs.BearingElement(2, kxx=stf, cxx=0) >>> rotor = rs.Rotor(shaft_elm, [disk0], [bearing0, bearing1]) >>> modal = rotor.run_modal(speed=0) >>> modal.wd[0] 215.3707...
Methods
- A(speed=0, frequency=None, synchronous=False)#
State space matrix for an instance of a rotor.
- Parameters:
- speed: float, optional
Rotor speed. Default is 0.
- frequencyfloat, optional
Excitation frequency. Default is rotor speed.
- synchronousbool, optional
If True a synchronous analysis is carried out. Default is False.
- Returns:
- Anp.ndarray
State space matrix for the rotor.
Examples
>>> rotor = rotor_example() >>> np.round(rotor.A()[75:83, :2]) + 0. array([[ 0., 10927.], [-10924., 0.], [ 0., 0.], [ -174., 0.], [ 0., -174.], [ 0., 0.], [ 0., 10723.], [-10719., 0.]])
- C(frequency)#
Damping matrix for an instance of a rotor.
- Parameters:
- frequencyfloat
Excitation frequency.
- Returns:
- C0np.ndarray
Damping matrix for the rotor.
Examples
>>> rotor = compressor_example() >>> rotor.C(0)[:4, :4] array([[0., 0., 0., 0.], [0., 0., 0., 0.], [0., 0., 0., 0.], [0., 0., 0., 0.]])
- G()#
Gyroscopic matrix for an instance of a rotor.
- Returns:
- G0np.ndarray
Gyroscopic matrix for the rotor.
Examples
>>> rotor = rotor_example() >>> rotor.G()[:4, :4] array([[ 0. , 0.01943344, 0. , -0.00022681], [-0.01943344, 0. , 0. , 0. ], [ 0. , 0. , 0. , 0. ], [ 0.00022681, 0. , 0. , 0. ]])
- K(frequency)#
Stiffness matrix for an instance of a rotor.
- Parameters:
- frequencyfloat, optional
Excitation frequency.
- Returns:
- K0np.ndarray
Stiffness matrix for the rotor.
Examples
>>> rotor = rotor_example() >>> np.round(rotor.K(0)[:4, :4] / 1e6) array([[ 4.700e+01, 0.000e+00, 0.000e+00, 0.000e+00], [ 0.000e+00, 4.600e+01, 0.000e+00, -6.000e+00], [ 0.000e+00, 0.000e+00, 1.657e+03, 0.000e+00], [ 0.000e+00, -6.000e+00, 0.000e+00, 1.000e+00]])
- Ksdt()#
Dynamic stiffness matrix for an instance of a rotor.
Stiffness matrix associated with the transient motion of the shaft and disks. It needs to be multiplied by the angular acceleration when considered in time dependent analyses.
- Returns:
- Ksdt0np.ndarray
Dynamic stiffness matrix for the rotor.
Examples
>>> rotor = rotor_example_6dof() >>> np.round(rotor.Ksdt()[:6, :6] * 1e3, 2) array([[ 0. , -23. , 0. , 0.48, 0. , 0. ], [ 0. , 0. , 0. , 0. , 0. , 0. ], [ 0. , 0. , 0. , 0. , 0. , 0. ], [ 0. , 0. , 0. , 0. , 0. , 0. ], [ 0. , -0.48, 0. , 0.16, 0. , 0. ], [ 0. , 0. , 0. , 0. , 0. , 0. ]])
- M(frequency=None, synchronous=False)#
Mass matrix for an instance of a rotor.
- Parameters:
- synchronousbool, optional
If True a synchronous analysis is carried out. Default is False.
- Returns:
- M0np.ndarray
Mass matrix for the rotor.
Examples
>>> rotor = rotor_example() >>> rotor.M(0)[:4, :4] array([[ 1.42050794, 0. , 0. , 0. ], [ 0. , 1.42050794, 0. , -0.04931719], [ 0. , 0. , 1.27790826, 0. ], [ 0. , -0.04931719, 0. , 0.00231392]])
- __init__(shaft_elements, disk_elements=None, bearing_elements=None, point_mass_elements=None, min_w=None, max_w=None, rated_w=None, tag=None)#
- add_nodes(new_nodes_pos)#
Add nodes to rotor.
This method returns the modified rotor with additional nodes according to the positions of the new nodes provided.
- Parameters:
- new_nodes_poslist
List with the position of the new nodes.
- Returns:
- A rotor object.
Examples
>>> import ross as rs >>> rotor = rs.rotor_example() >>> new_rotor = rotor.add_nodes([0.62, 1.11]) >>> shaft_elements = new_rotor.shaft_elements >>> len(shaft_elements) 8 >>> round(shaft_elements[3].L, 2) 0.13 >>> round(shaft_elements[6].L, 2) 0.14
- convergence(n_eigval=0, err_max=0.01)#
Run convergence analysis.
Function to analyze the eigenvalues convergence through the number of shaft elements. Every new run doubles the number os shaft elements.
- Parameters:
- n_eigvalint
The nth eigenvalue which the convergence analysis will run. Default is 0 (the first eigenvalue).
- err_maxfloat
Maximum allowable convergence error. Default is 1e-02
- Returns:
- resultsAn instance of ConvergenceResults class, which is used to post-process
- results. Attributes stored:
- el_numarray
Array with number of elements in each iteraction
- eigv_arrarray
Array with the n’th natural frequency in each iteraction
- error_arrarray
Array with the relative error in each iteraction
- classmethod from_section(leng_data, idl_data, odl_data, idr_data=None, odr_data=None, material_data=None, disk_data=None, brg_seal_data=None, min_w=None, max_w=None, rated_w=None, nel_r=1, tag=None)#
Build rotor from sections.
This class is an alternative to build rotors from separated sections. Each section has the same number (n) of shaft elements.
- Parameters:
- leng_datalist
List with the lengths of rotor regions.
- idl_datalist
List with the inner diameters of rotor regions (Left Station).
- odl_datalist
List with the outer diameters of rotor regions (Left Station).
- idr_datalist, optional
List with the inner diameters of rotor regions (Right Station). Default is equal to idl_data (cylindrical element).
- odr_datalist, optional
List with the outer diameters of rotor regions (Right Station). Default is equal to odl_data (cylindrical element).
- material_dataross.material or list of ross.material
Defines a single material for all sections or each section can be defined by a material individually.
- disk_datadict, optional
Dict holding disks datas. Example : disk_data=DiskElement.from_geometry(n=2,
material=steel, width=0.07, i_d=0, o_d=0.28 )
*See ‘disk_element.py’ docstring for more information*
- brg_seal_datadict, optional
Dict holding lists of bearings and seals datas. Example : brg_seal_data=BearingElement(n=1, kxx=1e6, cxx=0,
kyy=1e6, cyy=0, kxy=0, cxy=0, kyx=0, cyx=0)
*See ‘bearing_seal_element.py’ docstring for more information*
- nel_rint, optional
Number or elements per shaft region. Default is 1.
- tagstr
A tag for the rotor
- Returns:
- A rotor object
- Raises:
- ValueError
Error raised if lists size do not match.
- AttributeError
Error raised if the shaft material is not defined.
- gravitational_force(g=-9.8065, direction='y', M=None, num_dof=None)#
Compute the gravitational force vector for the system.
- Parameters:
- gfloat, optional
Acceleration due to gravity. Default is -9.8065 m/s².
- direction{“x”, “y”, “z”}, optional
Direction in which gravity acts. Default is “y”.
- Mndarray, optional
Mass matrix of the system. If None, the internal mass matrix is used.
- num_dofint, optional
Number of degrees of freedom per node. If None, the internal value is used.
- Returns:
- forcendarray
Gravitational force (weight) vector of shape (ndof,).
Examples
>>> rotor = compressor_example() >>> force = rotor.gravitational_force() >>> force[:4] array([ 0. , -3.12941854, 0. , 0.01851573])
- init_ambs_for_integrate(**kwargs)#
- integrate_system(speed, F, t, **kwargs)#
Time integration for a rotor system.
This method returns the time response for a rotor given a force, time and speed based on time integration with the Newmark method.
- Parameters:
- speedfloat or array_like
Rotor speed.
- Fndarray
Force array (needs to have the same length as time array).
- tndarray
Time array.
- **kwargsoptional
Additional keyword arguments can be passed to define the parameters of the Newmark method if it is used (e.g. gamma, beta, tol, …). See newmark for more details. Other optional arguments are listed below.
- model_reductiondict, optional
When model_reduction is provided, the corresponding reduction method is initialized. Dict keys:
- methodstr, optional
Reduction method to use, e.g., “guyan” or “pseudomodal”. Defaults to “guyan”.
- num_modesint, optional
Number of modes to reduce the model to, if pseudo-modal method is considered.
- include_nodeslist of int, optional
List of the nodes to be included, if Guyan reduction method is considered.
- dof_mappinglist of str, optional
List of the local DOFs to be considered when using Guyan reduction method. Valid values are: ‘x’, ‘y’, ‘z’, ‘alpha’, ‘beta’, ‘theta’, corresponding to:
‘x’ and ‘y’: lateral translations
‘z’: axial translation
‘alpha’: rotation about the x-axis
‘beta’: rotation about the y-axis
‘theta’: torsional rotation (about the z-axis)
Default is [‘x’, ‘y’].
- include_dofs (list of int, optional):
Additional degrees of freedom (DOFs) to include in the reduction, such as DOFs with applied forces or probe locations when using Guyan reduction method.
- add_to_RHScallable, optional
An optional function that computes and returns an additional array to be added to the right-hand side of the equation of motion. This function should take the time step number as argument, and can take optional arguments corresponding to the current state of the rotor system, including the displacements disp_resp, velocities velc_resp, and acceleration accl_resp. It should return an array of the same length as the degrees of freedom of the rotor system rotor.ndof. This function allows for the incorporation of supplementary terms or external effects in the rotor system dynamics beyond the specified force input during the time integration process.
- Returns:
- tndarray
Time values for the output.
- youtndarray
System response.
Examples
>>> import ross as rs >>> rotor = rs.compressor_example() >>> size = 10000 >>> node = 3 >>> speed = 500.0 >>> accel = 0.0 >>> t = np.linspace(0, 10, size) >>> F = np.zeros((size, rotor.ndof)) >>> F[:, rotor.number_dof * node + 0] = 10 * np.cos(2 * t) >>> F[:, rotor.number_dof * node + 1] = 10 * np.sin(2 * t) >>> t, yout = rotor.integrate_system(speed, F, t) Running direct method >>> yout[:, rotor.number_dof * node + 1] array([0.00000000e+00, 2.07239823e-10, 7.80952429e-10, ..., 1.21848307e-07, 1.21957287e-07, 1.22065778e-07])
- classmethod load(file)#
Load rotor from toml file.
- Parameters:
- filestr or pathlib.Path
String or Path for a .toml file.
- Returns:
- rotorross.rotor.Rotor
- magnetic_bearing_controller(step, magnetic_bearings, time_step, disp_resp, **kwargs)#
Compute control forces for Active Magnetic Bearings (AMBs).
This method calculates the magnetic control forces generated by active magnetic bearings (AMBs) at each time step using a PID control law. The forces are based on the measured displacements and can optionally include external disturbances for sensitivity analysis.
If sensitivity analysis is enabled via keyword arguments, the method injects a known disturbance at a specific DoF and logs excitation, disturbed, and sensor signals for post-processing.
- Parameters:
- stepint
Current time step index in the simulation.
- magnetic_bearingslist
List of MagneticBearingElement objects used for force computation.
- time_stepfloat
Time increment used in the numerical integration scheme (in seconds).
- disp_respndarray
Displacement response vector of the rotor at the current time step. The size must match the number of rotor DoFs.
- Returns:
- magnetic_forcendarray
Force vector containing control forces applied by each magnetic bearing in the rotor system. Has the same length as self.ndof.
- Other Parameters:
- sensitivity_compute_dofint, optional
Index of the DoF where a disturbance signal is applied (for sensitivity analysis).
- sensitivity_disturbancendarray, optional
Disturbance signal array (e.g., chirp) to be injected at the specified DoF.
- sensitivity_result_valuesdict, optional
- Dictionary to store the time history of:
“excitation_signal”
“disturbed_signal”
“sensor_signal”
for post-processing in sensitivity computations.
Notes
The control forces are applied in both x and y directions at each AMB location.
The actual PID computation is delegated to the compute_pid_amb function.
If sensitivity_compute_dof is provided, the excitation is applied to that DoF only.
Examples
>>> import ross as rs >>> import numpy as np >>> rotor = rs.rotor_amb_example() >>> dt, speed, step = 1e-4, 1000, 1 >>> t = np.arange(0, 5 * dt, dt) >>> node = [27, 29] >>> mass = [10, 10] >>> F = np.zeros((len(t), rotor.ndof)) >>> for n, m in zip(node,mass): ... F[:, 6 * n + 0] = m * np.cos((speed * t)) ... F[:, 6 * n + 1] = (m-5) * np.sin((speed * t)) >>> response = rotor.run_time_response(speed, F, t, method = "newmark") Running direct method >>> magnetic_bearings = [brg for brg in rotor.bearing_elements if isinstance(brg, rs.bearing_seal_element.MagneticBearingElement)] >>> magnetic_force = rotor.magnetic_bearing_controller(step, magnetic_bearings, dt, response.yout[-1,:]) >>> np.nonzero(magnetic_force)[0] array([ 72, 73, 258, 259]) >>> magnetic_force[np.nonzero(magnetic_force)[0]] array([-7.24276404e-04, -1.42153354e-05, -1.17641699e-04, 2.39844354e-05])
- plot_rotor(nodes=1, check_sld=False, length_units='m', **kwargs)#
Plot a rotor object.
This function will take a rotor object and plot its elements representation using Plotly.
- Parameters:
- nodesint, optional
Increment that will be used to plot nodes label.
- check_sldbool
If True, checks the slenderness ratio for each element. The shaft elements which has a slenderness ratio < 1.6 will be displayed in yellow color.
- length_unitsstr, optional
length units to length and diameter. Default is ‘m’.
- kwargsoptional
Additional key word arguments can be passed to change the plot layout only (e.g. width=1000, height=800, …). *See Plotly Python Figure Reference for more information.
- Returns:
- figplotly.graph_objects.Figure
The figure object with the rotor representation.
- run_amb_sensitivity(speed, t_max, dt, disturbance_amplitude=1e-05, disturbance_min_frequency=0.001, disturbance_max_frequency=150, amb_tags=None, sensors_theta=0.7853981633974483, verbose=1)#
Run Active Magnetic Bearing (AMB) sensitivity analysis.
This method performs a frequency-domain sensitivity analysis of the rotor system equipped with active magnetic bearings (AMBs). The analysis uses a logarithmic chirp excitation applied as an external disturbance force to compute the system’s frequency response at the AMB-controlled degrees of freedom (DoFs). The results provide magnitude and phase sensitivity functions for each AMB in both x and y directions.
- Parameters:
- speedfloat
Rotational speed of the rotor in rad/s.
- t_maxfloat
Total time duration of the simulation in seconds.
- dtfloat
Time step for the simulation in seconds.
- disturbance_amplitudefloat, optional
Amplitude of the excitation chirp signal applied as a disturbance. Default is 10e-6.
- disturbance_min_frequencyfloat, optional
Minimum frequency (in Hz) of the logarithmic chirp signal used for excitation. The chirp sweeps from this frequency up to disturbance_max_frequency. Default is 1e-3 Hz.
- disturbance_max_frequencyfloat, optional
Maximum frequency (in Hz) of the logarithmic chirp signal used for excitation. Default is 150 Hz.
- amb_tagslist of str, optional
List of magnetic bearing tags to include in the sensitivity analysis. If None or empty, all MagneticBearingElement instances in the rotor are used. If provided, only the AMBs matching the specified tags will be analyzed. Raises a RuntimeError if no AMB with the given tag is found.
- sensors_thetafloat, optional
Angular position of the Active Magnetic Bearing (AMB) sensors, in radians. This angle defines the orientation of the sensor coordinate system (v, w) relative to the global coordinate system (x, y). A positive angle corresponds to a counter-clockwise rotation. Default is 45 degrees (π/4 rad).
- verboseint, optional
Controls the verbosity of the method. If 1 or greater, both the simulation time and the forces produced by the AMBs are presented. If 0, no output is shown. Default is 1.
- Returns:
- resultsSensitivityResults
Object containing sensitivity magnitude, phase, and frequency vectors for each magnetic bearing tag and direction (‘x’, ‘y’). Also includes the excitation, disturbed, and sensor signals used in the computation.
Notes
The excitation is a logarithmic chirp sweeping from disturbance_min_frequency to disturbance_max_frequency (Hz).
The excitation is applied individually to each DoF controlled by an AMB.
The method assumes that the rotor contains MagneticBearingElement instances.
A Newmark time integration scheme is used internally via run_time_response().
Examples
>>> import ross as rs >>> rotor = rs.rotor_amb_example()
>>> # Run sensitivity for all magnetic bearings in the rotor (default sweep) >>> sensitivity_results = rotor.run_amb_sensitivity(speed=314.16, t_max=5e-4, dt=1e-4) Running direct method...
>>> # Run sensitivity only for a specific AMB tag (e.g., "Magnetic Bearing 0") >>> sensitivity_results = rotor.run_amb_sensitivity( ... speed=314.16, t_max=5e-4, dt=1e-4, amb_tags=["Magnetic Bearing 0"] ... ) Running direct method...
>>> # Run sensitivity with a custom chirp band (0.1 Hz to 200 Hz) >>> sensitivity_results = rotor.run_amb_sensitivity( ... speed=314.16, t_max=5e-4, dt=1e-4, ... disturbance_min_frequency=0.1, disturbance_max_frequency=200.0 ... ) Running direct method...
>>> # Accessing maximum absolute sensitivities for "Magnetic Bearing 0" >>> max_sens_bearing_0_x = sensitivity_results.max_abs_sensitivities["Magnetic Bearing 0"]["x"] >>> max_sens_bearing_0_y = sensitivity_results.max_abs_sensitivities["Magnetic Bearing 0"]["y"]
>>> # Plotting the sensitivities for all AMBs and axes >>> fig = sensitivity_results.plot( ... frequency_units="Hz", phase_unit="degree", ... magnitude_scale="decibel", xaxis_type="log" ... )
>>> # Plotting the time results used in sensitivity calculation >>> fig = sensitivity_results.plot_time_results()
- run_campbell(speed_range, frequencies=6, frequency_type='wd', torsional_analysis=False)#
Calculate the Campbell diagram.
This function will calculate the damped natural frequencies for a speed range.
- Available plotting methods:
.plot()
- Parameters:
- speed_rangearray, pint.Quantity
Array with the speed range in rad/s.
- frequenciesint, optional
Number of frequencies that will be calculated. Default is 6.
- frequency_typestr, optional
Choose between displaying results related to the undamped natural frequencies (“wn”) or damped natural frequencies (“wd”). The default is “wd”.
- torsional_analysisbool, optional
If True, performs a separate torsional analysis and returns the respective modes in the Campbell diagram. In this case, a system with only torsional degrees of freedom is considered, thus disregarding coupled modes (lateral + torsional). Default is False.
- Returns:
- resultsross.CampbellResults
For more information on attributes and methods available see:
ross.CampbellResults
Examples
>>> import ross as rs >>> rotor1 = rs.rotor_example() >>> speed = np.linspace(0, 400, 11)
Diagram with undamped natural frequencies >>> camp = rotor1.run_campbell(speed, frequency_type=”wn”)
Diagram with damped natural frequencies >>> camp = rotor1.run_campbell(speed)
Plotting Campbell Diagram >>> fig = camp.plot()
- run_crack(n, depth_ratio, node, unbalance_magnitude, unbalance_phase, speed, t, crack_model='Mayes', cross_divisions=None, **kwargs)#
Run analysis for the rotor system with crack given an unbalance force.
Crack object is instantiated and system’s time response is simulated.
- Parameters:
- nfloat
Element number where the crack is located.
- depth_ratiofloat
Crack depth ratio related to the diameter of the crack container element. A depth value of 0.1 is equal to 10%, 0.2 equal to 20%, and so on.
- nodelist, int
Node where the unbalance is applied.
- unbalance_magnitudelist, float, pint.Quantity
Unbalance magnitude (kg.m).
- unbalance_phaselist, float, pint.Quantity
Unbalance phase (rad).
- speedfloat or array_like, pint.Quantity
Rotor speed.
- Farray
Force array (needs to have the same number of rows as time array). Each column corresponds to a dof and each row to a time.
- tarray
Time array.
- crack_modelstring, optional
String containing type of crack model chosed. The available types are: “Mayes”, “Gasch”, “Flex Open” and “Flex Breathing”. Default is “Mayes”.
- cross_divisions: float, optional
Number of square divisions into which the cross-section of the cracked element will be divided in the analysis conducted for the Flex Breathing model.
- **kwargsoptional
Additional keyword arguments can be passed to define the parameters of the Newmark method if it is used (e.g. gamma, beta, tol, …). See ross.utils.newmark for more details. Other keyword arguments can also be passed to be used in numerical integration (e.g. model_reduction). See Rotor.integrate_system for more details.
- Returns:
- resultsross.TimeResponseResults
For more information on attributes and methods available see:
ross.TimeResponseResults
Examples
>>> import ross as rs >>> from ross.probe import Probe >>> from ross.units import Q_ >>> rotor = rs.rotor_example_with_damping() >>> n1 = rotor.disk_elements[0].n >>> n2 = rotor.disk_elements[1].n >>> results = rotor.run_crack( ... n=18, ... depth_ratio=0.2, ... node=[n1, n2], ... unbalance_magnitude=[5e-4, 0], ... unbalance_phase=[-np.pi / 2, 0], ... crack_model="Mayes", ... speed=Q_(1200, "RPM"), ... t=np.arange(0, 0.5, 0.0001), ... model_reduction={"num_modes": 12}, # Pseudo-modal method ... ) Running with model reduction: pseudomodal >>> probe1 = Probe(14, 0) >>> probe2 = Probe(22, 0) >>> fig1 = results.plot_1d([probe1, probe2]) >>> fig2 = results.plot_dfft( ... [probe1, probe2], ... frequency_range=Q_((0, 200), "Hz"), ... yaxis_type="log", ... )
- run_critical_speed(speed_range=None, num_modes=12, rtol=0.005)#
Calculate the critical speeds and damping ratios for the rotor model.
This function runs an iterative method over “run_modal()” to minimize (using scipy.optimize.newton) the error between the rotor speed and the rotor critical speeds (rotor speed - critical speed).
Differently from run_modal(), this function doesn’t take a speed input because it iterates over the natural frequencies calculated in the last iteration. The initial value is considered to be the undamped natural frequecies for speed = 0 (no gyroscopic effect).
Once the error is within an acceptable range defined by “rtol”, it returns the approximated critical speed.
With the critical speeds calculated, the function uses the results to calculate the log dec and damping ratios for each critical speed.
- Parameters:
- speed_rangetuple, optional, pint.Quantity
Tuple (start, end) with the desired range of frequencies (rad/s). The function returns all eigenvalues within this range.
- num_modesint, optional
The number of eigenvalues and eigenvectors to be calculated using ARPACK. If sparse=True, it determines the number of eigenvalues and eigenvectors to be calculated. It must be smaller than Rotor.ndof - 1. It is not possible to compute all eigenvectors of a matrix with ARPACK. If speed_range is not None, num_modes is overrided. Default is 12.
- rtolfloat, optional
Tolerance (relative) for termination. Applied to scipy.optimize.newton. Default is 0.005 (0.5%).
- Returns:
- resultsross.CriticalSpeedResults
For more information on attributes and methods available see:
ross.CriticalSpeedResults
Examples
>>> import ross as rs >>> rotor = rs.rotor_example()
Finding the first Nth critical speeds >>> results = rotor.run_critical_speed(num_modes=8) >>> np.round(results.wd()) array([ 92., 96., 271., 300.]) >>> np.round(results.wn()) array([ 92., 96., 271., 300.])
Finding the first critical speeds within a speed range >>> results = rotor.run_critical_speed(speed_range=(100, 1000)) >>> np.round(results.wd()) array([271., 300., 636., 774., 867.])
Changing output units >>> np.round(results.wd(“rpm”)) array([2590., 2868., 6074., 7394., 8278.])
Retrieving whirl directions >>> results.whirl_direction # doctest: +ELLIPSIS array([…
- run_forced_response(force=None, speed_range=None, modes=None, cluster_points=False, num_modes=12, num_points=10, rtol=0.005, unbalance=None)#
Forced response for a mdof system.
This method returns the unbalanced response for a mdof system given magnitude and phase of the unbalance, the node where it’s applied and a frequency range.
- Available plotting methods:
.plot() .plot_magnitude() .plot_phase() .plot_polar_bode() .plot_deflected_shape() .plot_bending_moment() .plot_deflected_shape_3d() .plot_deflected_shape_2d()
- Parameters:
- forcelist, array, pint.Quantity
Unbalance force in each degree of freedom for each value in omega
- speed_rangelist, array, pint.Quantity
Array with the desired range of frequencies
- modeslist, optional
Modes that will be used to calculate the frequency response (all modes will be used if a list is not given).
- unbalancearray, optional
Array with the unbalance data (node, magnitude and phase) to be plotted with deflected shape. This argument is set only if running an unbalance response analysis. Default is None.
- cluster_pointsbool, optional
boolean to activate the automatic frequency spacing method. If True, the method uses _clustering_points() to create an speed_range. Default is False
- num_pointsint, optional
The number of points generated per critical speed. The method set the same number of points for slightly less and slightly higher than the natural circular frequency. It means there’ll be num_points greater and num_points smaller than a given critical speed. num_points may be between 2 and 12. Anything above this range defaults to 10 and anything below this range defaults to 4. The default is 10.
- num_modes
The number of eigenvalues and eigenvectors to be calculated using ARPACK. It also defines the range for the output array, since the method generates points only for the critical speed calculated by run_critical_speed(). Default is 12.
- rtolfloat, optional
Tolerance (relative) for termination. Applied to scipy.optimize.newton to calculate the approximated critical speeds. Default is 0.005 (0.5%).
- Returns:
- resultsross.ForcedResponseResults
For more information on attributes and methods available see:
ross.ForcedResponseResults
Examples
>>> rotor = rotor_example() >>> speed = np.linspace(0, 1000, 101) >>> force = rotor._unbalance_force(3, 10.0, 0.0, speed) >>> resp = rotor.run_forced_response(force=force, speed_range=speed) >>> abs(resp.forced_resp) array([[0.00000000e+00, 5.06073311e-04, 2.10044826e-03, ...
Using clustered points option. Set cluster_points=True and choose how many modes the method must search and how many points to add just before and after each critical speed.
>>> response = rotor.run_forced_response( ... force=force, cluster_points=True, num_modes=12, num_points=5 ... ) >>> response.speed_range.shape (61,)
- run_freq_response(speed_range=None, modes=None, cluster_points=False, num_modes=12, num_points=10, rtol=0.005, free_free=False)#
Frequency response for a mdof system.
This method returns the frequency response for a mdof system given a range of frequencies and the modes that will be used.
- Available plotting methods:
.plot() .plot_magnitude() .plot_phase() .plot_polar_bode()
- Parameters:
- speed_rangearray, optional, pint.Quantity
Array with the desired range of frequencies. Default is 0 to 1.5 x highest damped natural frequency.
- modeslist, optional
Modes that will be used to calculate the frequency response (all modes will be used if a list is not given).
- cluster_pointsbool, optional
boolean to activate the automatic frequency spacing method. If True, the method uses _clustering_points() to create an speed_range. Default is False
- num_pointsint, optional
The number of points generated per critical speed. The method set the same number of points for slightly less and slightly higher than the natural circular frequency. It means there’ll be num_points greater and num_points smaller than a given critical speed. num_points may be between 2 and 12. Anything above this range defaults to 10 and anything below this range defaults to 4. The default is 10.
- num_modes
The number of eigenvalues and eigenvectors to be calculated using ARPACK. It also defines the range for the output array, since the method generates points only for the critical speed calculated by run_critical_speed(). Default is 12.
- rtolfloat, optional
Tolerance (relative) for termination. Applied to scipy.optimize.newton to calculate the approximated critical speeds. Default is 0.005 (0.5%).
- free_freebool, optional
If True, the method will consider the rotor system as free-free. Default is False.
- Returns:
- resultsross.FrequencyResponseResults
For more information on attributes and methods available see:
ross.FrequencyResponseResults
Examples
>>> import ross as rs >>> rotor = rs.rotor_example() >>> speed =np.linspace(0, 1000, 101) >>> response = rotor.run_freq_response(speed_range=speed)
Return the response amplitude >>> abs(response.freq_resp) # doctest: +ELLIPSIS array([[[0.00000000e+00, 1.00261725e-06, 1.01076952e-06, …
Return the response phase >>> np.angle(response.freq_resp) # doctest: +ELLIPSIS array([[[…
Using clustered points option. Set cluster_points=True and choose how many modes the method must search and how many points to add just before and after each critical speed.
>>> response = rotor.run_freq_response(cluster_points=True, num_points=5) >>> response.speed_range.shape (61,)
Selecting the disirable modes, if you want a reduced model: >>> response = rotor.run_freq_response(speed_range=speed, modes=[0, 1, 2, 3, 4]) >>> abs(response.freq_resp) # doctest: +ELLIPSIS array([[[0.00000000e+00, 1.00261725e-06, 1.01076952e-06, …
Plotting frequency response function: >>> fig = response.plot(inp=13, out=13)
To plot velocity and acceleration responses, you must change amplitude_units from “[length]/[force]” units to “[speed]/[force]” or “[acceleration]/[force]” respectively
Plotting velocity response >>> fig = response.plot(inp=13, out=13, amplitude_units=”m/s/N”)
Plotting acceleration response >>> fig = response.plot(inp=13, out=13, amplitude_units=”m/s**2/N”)
- run_level1(n=5, stiffness_range=None, num=5, **kwargs)#
Plot level 1 stability analysis.
This method will plot the stability 1 analysis for a given stiffness range.
- Parameters:
- nint
Number of steps in the range. Default is 5.
- stiffness_rangetuple, optional
Tuple with (start, end) for stiffness range. This will be used to create an evenly numbers spaced evenly on a log scale to create a better visualization (see np.logspace).
- kwargsoptional
Additional key word arguments can be passed to change the plot layout only (e.g. width=1000, height=800, …). *See Plotly Python Figure Reference for more information.
- Returns:
- resultsross.Level1Results
For more information on attributes and methods available see:
ross.Level1Results
- run_misalignment(node, unbalance_magnitude, unbalance_phase, speed, t, coupling='flex', **kwargs)#
Run analysis for the rotor system with misalignment given an unbalance force.
Misalignment object is instantiated and system’s time response is simulated. There are two types of coupling: flexible (flex) and rigid, each with distinct parameters. These parameters are passed to the respective method through **kwargs.
- Parameters:
- nodelist, int
Node where the unbalance is applied.
- unbalance_magnitudelist, float, pint.Quantity
Unbalance magnitude (kg.m).
- unbalance_phaselist, float, pint.Quantity
Unbalance phase (rad).
- speedfloat or array_like, pint.Quantity
Rotor speed.
- Farray
Force array (needs to have the same number of rows as time array). Each column corresponds to a dof and each row to a time.
- tarray
Time array.
- couplingstr
Coupling type. The avaible types are: “flex” and “rigid”. Default is “flex”.
- **kwargsdictionary
- If coupling = “flex”, **kwargs receives:
- nfloat
Number of shaft element where the misalignment is ocurring.
- mis_type: string
Name of the chosen misalignment type. The avaible types are: “parallel”, “angular” and “combined”.
- mis_distance_xfloat, pint.Quantity
Parallel misalignment distance between driving rotor and driven rotor along X direction.
- mis_distance_yfloat, pint.Quantity
Parallel misalignment distance between driving rotor and driven rotor along Y direction.
- mis_anglefloat, pint.Quantity
Angular misalignment angle.
- radial_stiffnessfloat, pint.Quantity
Radial stiffness of flexible coupling.
- bending_stiffnessfloat, pint.Quantity
Bending stiffness of flexible coupling. Provide if mis_type is “angular” or “combined”.
- input_torquefloat, pint.Quantity
Driving torque. Default is 0.
- load_torquefloat, pint.Quantity
Driven torque. Default is 0.
- If coupling = “rigid”, **kwargs receives:
- nfloat
Number of shaft element where the misalignment is ocurring.
- mis_distancefloat, pint.Quantity
Parallel misalignment distance between driving rotor and driven rotor.
- input_torquefloat, pint.Quantity
Driving torque. Default is 0.
- load_torquefloat, pint.Quantity
Driven torque. Default is 0.
Additional keyword arguments can be passed to define the parameters of the Newmark method if it is used (e.g. gamma, beta, tol, …). See ross.utils.newmark for more details. Other keyword arguments can also be passed to be used in numerical integration (e.g. model_reduction). See Rotor.integrate_system for more details.
- Returns:
- resultsross.TimeResponseResults
For more information on attributes and methods available see:
ross.TimeResponseResults
Examples
>>> import ross as rs >>> from ross.probe import Probe >>> from ross.units import Q_ >>> rotor = rotor_example_with_damping() >>> n1 = rotor.disk_elements[0].n >>> n2 = rotor.disk_elements[1].n >>> results = rotor.run_misalignment( ... node=[n1, n2], ... unbalance_magnitude=[5e-4, 0], ... unbalance_phase=[-np.pi / 2, 0], ... speed=Q_(1200, "RPM"), ... t=np.arange(0, 0.5, 0.0001), ... coupling="rigid", ... n=0, ... mis_distance=2e-4, ... input_torque=0, ... load_torque=0, ... model_reduction={"num_modes": 12}, # Pseudo-modal method ... ) Running with model reduction: pseudomodal >>> probe1 = Probe(14, 0) >>> probe2 = Probe(22, 0) >>> fig1 = results.plot_1d([probe1, probe2]) >>> fig2 = results.plot_dfft( ... [probe1, probe2], ... frequency_range=Q_((0, 200), "Hz"), ... yaxis_type="log", ... )
- run_rubbing(n, distance, contact_stiffness, contact_damping, friction_coeff, node, unbalance_magnitude, unbalance_phase, speed, t, torque=False, **kwargs)#
Run analysis for the rotor system with rubbing given an unbalance force.
Rubbing object is instantiated and system’s time response is simulated.
- Parameters:
- nint
Number of shaft element where rubbing is ocurring.
- distancefloat, pint.Quantity
Distance between the housing and shaft surface.
- contact_stiffnessfloat, pint.Quantity
Contact stiffness.
- contact_dampingfloat, pint.Quantity
Contact damping.
- friction_coefffloat
Friction coefficient.
- nodelist, int
Node where the unbalance is applied.
- unbalance_magnitudelist, float, pint.Quantity
Unbalance magnitude (kg.m).
- unbalance_phaselist, float, pint.Quantity
Unbalance phase (rad).
- speedfloat or array_like, pint.Quantity
Rotor speed.
- Farray
Force array (needs to have the same number of rows as time array). Each column corresponds to a dof and each row to a time.
- tarray
Time array.
- torquebool, optional
If True a torque is considered by rubbing. Default is False.
- **kwargsoptional
Additional keyword arguments can be passed to define the parameters of the Newmark method if it is used (e.g. gamma, beta, tol, …). See ross.utils.newmark for more details. Other keyword arguments can also be passed to be used in numerical integration (e.g. model_reduction). See Rotor.integrate_system for more details.
- Returns:
- resultsross.TimeResponseResults
For more information on attributes and methods available see:
ross.TimeResponseResults
Examples
>>> import ross as rs >>> from ross.units import Q_ >>> from ross.probe import Probe >>> rotor = rotor_example_with_damping() >>> n1 = rotor.disk_elements[0].n >>> n2 = rotor.disk_elements[1].n >>> results = rotor.run_rubbing( ... n=12, ... distance=7.95e-5, ... contact_stiffness=1.1e6, ... contact_damping=40, ... friction_coeff=0.3, ... torque=False, ... node=[n1, n2], ... unbalance_magnitude=[5e-4, 0], ... unbalance_phase=[-np.pi / 2, 0], ... speed=Q_(1200, "RPM"), ... t=np.arange(0, 0.5, 0.0001), ... model_reduction={"num_modes": 12}, # Pseudo-modal method ... ) Running with model reduction: pseudomodal >>> probe1 = Probe(14, 0) >>> probe2 = Probe(22, 0) >>> fig1 = results.plot_1d([probe1, probe2]) >>> fig2 = results.plot_dfft( ... [probe1, probe2], ... frequency_range=Q_((0, 200), "Hz"), ... yaxis_type="log", ... )
- run_static()#
Run static analysis.
Static analysis calculates free-body diagram, deformed shaft, shearing force diagram and bending moment diagram.
- Available plotting methods:
.plot_deformation() .plot_bending_moment() .plot_shearing_force() .plot_free_body_diagram()
- Attributes:
- shaft_weight: float
Shaft total weight
- disk_forces_nodaldict
Relates the static force at each node due to the weight of disks
- bearing_forces_nodaldict
Relates the static force at each node due to the bearing reaction forces.
- bearing_forces_tagdict
Indicates the reaction force exerted by each bearing.
- disk_forces_tagdict
Indicates the force exerted by each disk.
- displacement_y: array
The shaft static displacement vector,
- Vx: array
Shearing force vector
- Bm: array
Bending moment vector
- Returns:
- resultsross.StaticResults
For more information on attributes and methods available see:
ross.StaticResults
- Raises:
- ValueError
Error raised if the rotor has no bearing elements.
- run_time_response(speed, F, t, method='default', **kwargs)#
Calculate the time response.
This function will take a rotor object and calculate its time response given a force and a time.
- Available plotting methods:
.plot_1d() .plot_2d() .plot_3d()
- Parameters:
- speedfloat or array_like, pint.Quantity
Rotor speed. Automatically, the Newmark method is chosen if speed has an array_like type.
- Farray
Force array (needs to have the same number of rows as time array). Each column corresponds to a dof and each row to a time.
- tarray
Time array.
- methodstr, optional
The Newmark method can be chosen by setting method=’newmark’.
- **kwargsoptional
Additional keyword arguments can be passed to define the parameters of the Newmark method if it is used (e.g. gamma, beta, tol, …). See ross.utils.newmark for more details. Other keyword arguments can also be passed to be used in numerical integration (e.g. model_reduction, add_to_RHS). See Rotor.integrate_system for more details.
- Returns:
- resultsross.TimeResponseResults
For more information on attributes and methods available see:
ross.TimeResponseResults
Examples
>>> from ross.probe import Probe >>> rotor = rotor_example() >>> speed = 500.0 >>> size = 1000 >>> node = 3 >>> probe1 = Probe(3, 0) >>> t = np.linspace(0, 10, size) >>> F = np.zeros((size, rotor.ndof)) >>> F[:, rotor.number_dof * node + 0] = 10 * np.cos(2 * t) >>> F[:, rotor.number_dof * node + 1] = 10 * np.sin(2 * t) >>> response = rotor.run_time_response(speed, F, t) >>> response.yout[:, rotor.number_dof * node + 1] array([ 0.00000000e+00, 1.86686693e-07, 8.39130663e-07, ... >>> # plot time response for a given probe: >>> fig1 = response.plot_1d(probe=[probe1]) >>> # plot orbit response - plotting 2D nodal orbit: >>> fig2 = response.plot_2d(node=node) >>> # plot orbit response - plotting 3D orbits - full rotor model: >>> fig3 = response.plot_3d()
- run_ucs(stiffness_range=None, bearing_frequency_range=None, num_modes=16, num=20, synchronous=False, **kwargs)#
Run Undamped Critical Speeds analyzes.
This method will run the undamped critical speed analyzes for a given range of stiffness values. If the range is not provided, the bearing stiffness at rated speed will be used to create a range.
- Parameters:
- stiffness_rangetuple, optional
Tuple with (start, end) for stiffness range in a log scale. In linear space, the sequence starts at
base ** start(base to the power of start) and ends withbase ** stop(see endpoint below). Here base is 10.0.- bearing_frequency_rangetuple, optional
The bearing frequency range used to calculate the intersection points. In some cases bearing coefficients will have to be extrapolated. The default is None. In this case the bearing frequency attribute is used.
- num_modesint, optional
Number of modes to be calculated. This uses scipy.sparse.eigs method. Default is 16. In this case 4 modes are plotted, since for each pair of eigenvalues calculated we have one wn, and we show only the forward mode in the plots.
- numint
Number of steps in the range. Default is 20.
- synchronousbool, optional
If True a synchronous analysis is carried out according to []. Default is False.
- Returns:
- resultsross.UCSResults
For more information on attributes and methods available see:
ross.UCSResults
- run_unbalance_response(node, unbalance_magnitude, unbalance_phase, frequency=None, modes=None, cluster_points=False, num_modes=12, num_points=10, rtol=0.005)#
Unbalanced response for a mdof system.
This method returns the unbalanced response for a mdof system given magnitide and phase of the unbalance, the node where it’s applied and a frequency range.
- Available plotting methods:
.plot() .plot_magnitude() .plot_phase() .plot_polar_bode() .plot_deflected_shape() .plot_bending_moment() .plot_deflected_shape_3d() .plot_deflected_shape_2d()
- Parameters:
- nodelist, int
Node where the unbalance is applied.
- unbalance_magnitudelist, float, pint.Quantity
Unbalance magnitude (kg.m).
- unbalance_phaselist, float, pint.Quantity
Unbalance phase (rad).
- frequencylist, pint.Quantity
List with the desired range of frequencies (rad/s).
- modeslist, optional
Modes that will be used to calculate the frequency response (all modes will be used if a list is not given).
- cluster_pointsbool, optional
boolean to activate the automatic frequency spacing method. If True, the method uses _clustering_points() to create an speed_range. Default is False
- num_pointsint, optional
The number of points generated per critical speed. The method set the same number of points for slightly less and slightly higher than the natural circular frequency. It means there’ll be num_points greater and num_points smaller than a given critical speed. num_points may be between 2 and 12. Anything above this range defaults to 10 and anything below this range defaults to 4. The default is 10.
- num_modes
The number of eigenvalues and eigenvectors to be calculated using ARPACK. It also defines the range for the output array, since the method generates points only for the critical speed calculated by run_critical_speed(). Default is 12.
- rtolfloat, optional
Tolerance (relative) for termination. Applied to scipy.optimize.newton to calculate the approximated critical speeds. Default is 0.005 (0.5%).
- Returns:
- resultsross.ForcedResponseResults
For more information on attributes and methods available see:
ross.ForcedResponseResults
Examples
>>> import ross as rs >>> rotor = rs.rotor_example() >>> speed = np.linspace(0, 1000, 101) >>> response = rotor.run_unbalance_response(node=3, ... unbalance_magnitude=10.0, ... unbalance_phase=0.0, ... frequency=speed)
Return the response amplitude >>> abs(response.forced_resp) # doctest: +ELLIPSIS array([[0.00000000e+00, 5.06073311e-04, 2.10044826e-03, …
Return the response phase >>> np.angle(response.forced_resp) # doctest: +ELLIPSIS array([[ 0. , 0. , 0. , …
Using clustered points option. Set cluster_points=True and choose how many modes the method must search and how many points to add just before and after each critical speed.
>>> response2 = rotor.run_unbalance_response( ... node=3, unbalance_magnitude=0.01, unbalance_phase=0.0, cluster_points=True, num_points=5 ... ) >>> response2.speed_range.shape (61,)
plot unbalance response: >>> probe_node = 3 >>> probe_angle = np.pi / 2 >>> probe_tag = “my_probe” # optional >>> fig = response.plot(probe=[rs.Probe(probe_node, probe_angle, tag=probe_tag)])
plot response for major or minor axis: >>> probe_node = 3 >>> probe_angle = “major” # for major axis >>> # probe_angle = “minor” # for minor axis >>> probe_tag = “my_probe” # optional >>> fig = response.plot(probe=[rs.Probe(probe_node, probe_angle, tag=probe_tag)])
To plot velocity and acceleration responses, you must change amplitude_units from “[length]” units to “[length]/[time]” or “[length]/[time] ** 2” respectively Plotting velocity response >>> fig = response.plot( … probe=[rs.Probe(probe_node, probe_angle)], … amplitude_units=”m/s” … )
Plotting acceleration response >>> fig = response.plot( … probe=[rs.Probe(probe_node, probe_angle)], … amplitude_units=”m/s**2” … )
Plotting deflected shape configuration Speed value must be in speed_range. >>> value = 600 >>> fig = response.plot_deflected_shape(speed=value)
- save(file)#
Save the rotor to a .toml file.
- Parameters:
- filestr or pathlib.Path
Examples
>>> from tempfile import tempdir >>> from pathlib import Path >>> # create path for temporary file >>> file = Path(tempdir) / 'rotor.toml' >>> rotor = rotor_example() >>> rotor.save(file)
- save_mat(file, speed, frequency=None)#
Save matrices and rotor model to a .mat file.
- Parameters:
- filestr, pathlib.Path
- speed: float
Rotor speed.
- frequency: float, optional
Excitation frequency. Default is rotor speed.
Examples
>>> from tempfile import tempdir >>> from pathlib import Path >>> # create path for temporary file >>> file = Path(tempdir) / 'new_matrices' >>> rotor = rotor_example() >>> rotor.save_mat(file, speed=0)
- summary()#
Plot the rotor summary.
This functioncreates a summary of the main parameters and attributes of the rotor model. The data is presented in a table format.
- Returns:
- resultsross.SummaryResults class
An instance of SumarryResults class to build the summary table
Examples
>>> rotor = rotor_example() >>> table = rotor.summary().plot() >>> # to display the plot use the command: >>> # show(table)
- time_response(speed, F, t, ic=None, method='default', **kwargs)#
Time response for a rotor.
This method returns the time response for a rotor given a force, time and initial conditions.
- Parameters:
- speedfloat or array_like
Rotor speed. Automatically, the Newmark method is chosen if speed has an array_like type.
- Farray
Force array (needs to have the same length as time array).
- tarray
Time array. (must have the same length than lti.B matrix)
- icarray, optional
The initial conditions on the state vector (zero by default).
- methodstr, optional
The Newmark method can be chosen by setting method=’newmark’.
- **kwargsoptional
Additional keyword arguments can be passed to define the parameters of the Newmark method if it is used (e.g. gamma, beta, tol, …). See ross.utils.newmark for more details. Other keyword arguments can also be passed to be used in numerical integration (e.g. model_reduction, add_to_RHS). See Rotor.integrate_system for more details.
- Returns:
- tarray
Time values for the output.
- youtarray
System response.
- xoutarray
Time evolution of the state vector.
Examples
>>> rotor = rotor_example() >>> speed = 0 >>> size = 28 >>> t = np.linspace(0, 5, size) >>> F = np.ones((size, rotor.ndof)) >>> rotor.time_response(speed, F, t) (array([0. , 0.18518519, 0.37037037, ...
- classmethod to_ross_only(rotor)#
Convert rotor with rsxl objects to ross only.
- transfer_matrix(speed=None, frequency=None, modes=None)#
Calculate the fer matrix for the frequency response function (FRF).
- Returns:
- Hmatrix
System transfer matrix
- unbalance_force_over_time(node, magnitude, phase, omega, t, return_all=False)#
Calculate unbalance forces for each time step.
This auxiliary function calculates the unbalanced forces by taking into account the magnitude and phase of the force. It generates an array of force values at each degree of freedom for the specified nodes at each time step, while also considering a range of frequencies.
- Parameters:
- nodelist, int
Nodes where the unbalance is applied.
- magnitudelist, float
Unbalance magnitude (kg.m) for each node.
- phaselist, float
Unbalance phase (rad) for each node.
- omegafloat, np.darray
Constant velocity or desired range of velocities (rad/s).
- tnp.darray
Time array (s).
- return_allbool, optional
If True, returns F0, theta, omega, and alpha. If False, returns only F0. Default is False.
- Returns:
- F0np.ndarray
Unbalance force at each degree of freedom for each time step.
- thetanp.ndarray
Angular positions for each time step.
- omeganp.ndarray
Angular velocities for each time step.
- alphanp.ndarray
Angular accelerations for each time step.
Examples
>>> rotor = rotor_example() >>> t = np.linspace(0, 10, 31) >>> omega = np.linspace(0, 1000, 31) >>> F = rotor.unbalance_force_over_time([3], [10.0], [0.0], omega, t) >>> F[18, :3] array([ 0. , 7632.15353293, -43492.18127561])
Attributes