LinearBeam
- class sharpy.linear.assembler.linearbeam.LinearBeam[source]
State space member
Define class for linear state-space realisation of GEBM flexible-body equations from SHARPy
timestep_info
class and with the nonlinear structural information.State-space models can be defined in continuous or discrete time (dt required). Modal projection, either on the damped or undamped modal shapes, is also available.
The beam state space has information on the states which will depend on whether the system is modal or expressed in physical coordinates.
If
modal
the state variables will beq
andq_dot
representing the modal displacements and the time derivatives.If
nodal
and free-flying, the state variables will beeta
for the flexible degrees of freedom (displacements and CRVs for each node (dim6)),V
representing the linear velocities at the A frame (dim3),W
representing the angular velocities at theA
frame (dim3), andorient
representing the orientation variable of theA
frame with respect toG
Notes on the settings:
modal_projection={True,False}
: determines whether to project the statesonto modal coordinates. Projection over damped or undamped modal shapes can be obtained selecting:
proj_modes = {'damped','undamped'}
while
inout_coords={'modes','nodal'}
determines whether the modal state-space inputs/outputs are modal coords or nodal degrees-of-freedom. If
modes
is selected, theKin
andKout
gain matrices are generated to transform nodal to modal dofs
dlti={True,False}
: if true, generates discrete-time system.The continuous to discrete transformation method is determined by:
discr_method={ 'newmark', # Newmark-beta 'zoh', # Zero-order hold 'bilinear'} # Bilinear (Tustin) transformation
DLTIs can be obtained directly using the Newmark-\(\beta\) method
discr_method='newmark'
newmark_damp=xx
withxx<<1.0
for full-states descriptions (
modal_projection=False
) and modal projection over the undamped structural modes (modal_projection=True
andproj_modes
). The Zero-order holder and bilinear methods, instead, work in all descriptions, but require the continuous state-space equations.
The settings that this solver accepts are given by a dictionary, with the following key-value pairs:
Name
Type
Description
Default
Options
modal_projection
bool
Use modal projection
True
inout_coords
str
Beam state space input/output coordinates
nodes
nodes
,modes
num_modes
int
Number of modes to retain
10
discrete_time
bool
Assemble beam in discrete time
True
dt
float
Discrete time system integration time step
0.001
proj_modes
str
Use
undamped
ordamped
modesundamped
damped
,undamped
discr_method
str
Discrete time assembly system method:
newmark
newmark
,zoh
,bilinear
newmark_damp
float
Newmark damping value. For systems assembled using
newmark
0.0001
use_euler
bool
Use euler angles for rigid body parametrisation
True
print_info
bool
Display information on screen
True
gravity
bool
Linearise gravitational forces
False
remove_dofs
list(str)
Remove desired degrees of freedom (flexible DOFs, linear velocities, rotational velocities, orientation)
[]
eta
,V
,W
,orient
remove_sym_modes
bool
Remove symmetric modes if wing is clamped
False
remove_rigid_states
bool
(For Stability Derivatives) - Remove RIGID STATES from SS leaving input/output channels unchanged
False
- assemble(t_ref=None)[source]
Assemble the beam state-space system.
- Parameters:
t_ref (float) – Scaling factor to non-dimensionalise the beam’s time step.
Returns:
- recover_accelerations(full_ss)[source]
For a system with displacement and velocity outputs (
full_ss
), recover the accelerations and append them as new output channels.This function produces an output gain that should then be connected in series to the desired system
- Parameters:
full_ss (libss.StateSpace) – State space for which to provide output gain to recover accelerations
- Returns:
Gain adding the accelerations as new output channels
- Return type:
libss.Gain
- remove_symmetric_modes()[source]
Removes symmetric modes when the wing is clamped at the midpoint.
It will force the wing tip displacements in
z
to be postive for all modes.Updates the mode shapes matrix, the natural frequencies and the number of modes.
- trim_nodes(trim_list=<class 'list'>)[source]
Removes degrees of freedom from the second order system.
- Parameters:
trim_list (list) – List of degrees of freedom to remove
eta
,V
,W
ororient
- unpack_flex_dof(eta, eta_dot=None)[source]
Unpacks a vector of structural displacements and velocities into a SHARPy familiar form of pos, psi and their time derivatives
- Parameters:
eta (np.array) – Vector of structural displacements
eta_dot (np.array (Optional) – Vector of structural velocities
- Returns:
- Containing
pos
,psi
,pos_dot
,psi_dot
ifeta_dot
is provided, else only the displacements are returned
- Containing
- Return type:
tuple
- unpack_ss_vector(x_n, u_n, struct_tstep)[source]
Warning
- Under development. Missing:
Accelerations
Double check the cartesian rotation vector
Tangential operator for the moments
Takes the state \(x = [\eta, \dot{\eta}]\) and input vectors \(u = N\) of a linearised beam and returns a SHARPy timestep instance, including the reference values.
- Parameters:
x_n (np.ndarray) – Structural beam state vector in nodal space
y_n (np.ndarray) – Beam input vector (nodal forces)
struct_tstep (utils.datastructures.StructTimeStepInfo) – Reference timestep used for linearisation
- Returns:
new timestep with linearised values added to the reference value
- Return type: