FlexDynamic

class sharpy.linear.src.lingebm.FlexDynamic(tsinfo, structure=None, custom_settings={})[source]

Define class for linear state-space realisation of GEBM flexible-body equations from SHARPy``timestep_info`` class and with the nonlinear structural information.

The linearised beam takes the following arguments:

Parameters:
  • tsinfo (sharpy.utils.datastructures.StructImeStepInfo) – Structural timestep containing the modal information

  • structure (sharpy.solvers.beamloader.Beam) – Beam class with the structural information

  • custom_settings (dict) – settings for the linearised beam

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 avaiable. The rad/s array wv can be optionally passed for freq. response analysis

To produce the state-space equations:

  1. Set the settings:
    1. modal_projection={True,False}: determines whether to project the states

      onto 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, the Kin and Kout gain matrices are generated to transform nodal to modal dofs

    2. 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 with xx<<1.0

      for full-states descriptions (modal_projection=False) and modal projection over the undamped structural modes (modal_projection=True and proj_modes). The Zero-order holder and bilinear methods, instead, work in all descriptions, but require the continuous state-space equations.

  2. Generate an instance of the beam

2. Run self.assemble(). The method accepts an additional parameter, Nmodes, which allows using a lower number of modes than specified in self.Nmodes

Examples

>>> beam_settings = {'modal_projection': True,
>>>             'inout_coords': 'modes',
>>>             'discrete_time': False,
>>>             'proj_modes': 'undamped',
>>>             'use_euler': True}
>>>
>>> beam = lingebm.FlexDynamic(tsstruct0, structure=data.structure, custom_settings=beam_settings)
>>>
>>> beam.assemble()

Notes

  • Modal projection will automatically select between damped/undamped modes shapes, based on the data available from tsinfo.

  • If the full system matrices are available, use the modal_sol methods to override mode-shapes and eigenvectors

assemble(Nmodes=None)[source]

Assemble state-space model

Several assembly options are available:

  1. Discrete-time, Newmark-\(\beta\):
    • Modal projection onto undamped modes. It uses the modal projection such that the generalised coordinates \(\eta\) are transformed into modal space by

      \[\mathbf{\eta} = \mathbf{\Phi\,q}\]

      where \(\mathbf{\Phi}\) are the first Nmodes right eigenvectors. Therefore, the equation of motion can be re-written such that the modes normalise the mass matrix to become the identity matrix.

      \[\mathbf{I_{Nmodes}}\mathbf{\ddot{q}} + \mathbf{\Lambda_{Nmodes}\,q} = 0\]

      The system is then assembled in Newmark-\(\beta\) form as detailed in newmark_ss()

    • Full size system assembly. No modifications are made to the mass, damping or stiffness matrices and the system is directly assembled by newmark_ss().

  2. Continuous time state-space

Parameters:

Nmodes (int) – number of modes to retain

cont2disc(dt=None)[source]

Convert continuous-time SS model into

converge_modal(wv=None, tol=None, Yref=None, Print=False)[source]

Determine number of modes required to achieve a certain convergence of the modal solution in a prescribed frequency range wv. The H-infinity norm of the error w.r.t. Yref is used for assessing convergence.

Warning

if a reference freq. response, Yref, is not provided, the full- state continuous-time frequency response is used as reference. This requires the full-states matrices Mstr, Cstr, Kstr to be available.

euler_propagation_equations(tsstr)[source]

Introduce the linearised Euler propagation equations that relate the body fixed angular velocities to the Earth fixed Euler angles.

This method will remove the quaternion propagation equations created by SHARPy; the resulting system will have 9 rigid degrees of freedom.

Parameters:

tsstr

Returns:

freqresp(wv=None, bode=True)[source]

Computes the frequency response of the current state-space model. If self.modal=True, the in/out are determined according to self.inout_coords

linearise_applied_forces(tsstr=None)[source]

Linearise externally applied follower forces given in the local B reference frame.

Updates the stiffness matrix with terms arising from this linearisation.

The linearised beam equations are expressed in the following frames of reference:

  • Nodal forces: \(\delta \mathbf{f}_A\)

  • Nodal moments: \(\delta(T^T \mathbf{m}_B)\)

  • Total forces (rigid body equations): \(\delta \mathbf{F}_A\)

  • Total moments (rigid body equations): \(\delta \mathbf{M}_A\)

Thus, when linearising externally applied follower forces projected onto the appropriate frame

\[\boldsymbol{f}_A^{ext} = C^{AB}(\boldsymbol{\psi})\boldsymbol{f}^{ext}_B\]

the following terms appear:

\[\delta\boldsymbol{f}_A^{ext} = \frac{\partial}{\partial\boldsymbol{\psi}} \left(C^{AB}(\boldsymbol{\psi})\boldsymbol{f}^{ext}_{0,B}\right)\delta\boldsymbol{\psi} + C^{AB}_0\delta\boldsymbol f_B^{ext}\]

where the \(\delta\boldsymbol{\psi}\) is a stiffenning term that needs to be included in the stiffness matrix. The terms will appear in the rows relating to the translational degrees of freedom and the columns that correspond to the cartesian rotation vector.

\[K_{ss}^{f,\Psi} \leftarrow -\frac{\partial}{\partial\boldsymbol{\psi}} \left(C^{AB}(\boldsymbol{\psi})\boldsymbol{f}^{ext}_{0,B}\right)\]

Externally applied moments in the material frame \(\boldsymbol{m}_B^{ext}\) result in the following linearised expression:

\[\delta(T^\top\boldsymbol{m}_B) = \frac{\partial}{\partial\boldsymbol{\psi}}\left( T^\top(\boldsymbol{\psi})\boldsymbol{m}^{ext}_{0,B}\right)\delta\boldsymbol{\psi} + T_0^\top \delta\boldsymbol{m}_B^{ext}\]

Which results in the following stiffenning term:

\[K_{ss}^{m,\Psi} \leftarrow -\frac{\partial}{\partial\boldsymbol{\psi}}\left( T^\top(\boldsymbol{\psi})\boldsymbol{m}^{ext}_{0,B}\right)\]

The total contribution of moments must be summed up for the rigid body equations, and include contributions due to externally applied forces as well as moments:

\[\boldsymbol{M}_A^{ext} = \sum_n \tilde{\boldsymbol{R}}_A C^{AB}(\boldsymbol{\psi}) \boldsymbol{f}_B^{ext} + \sum C^{AB}(\boldsymbol{\psi})\boldsymbol{m}_B^{ext}\]

The linearisation of this term becomes

\[\delta\boldsymbol{M}_A^{ext} = \sum\left(-\widetilde{C^{AB}_0 \boldsymbol{f}_{0,B}^{ext}}\delta \boldsymbol{R}_A + \widetilde{\boldsymbol{R}}\frac{\partial}{\partial\boldsymbol{\psi}}\left(C^{AB}\boldsymbol{f}_B\right) \delta \boldsymbol{\psi} + \widetilde{\boldsymbol{R}}C^{AB}\delta\boldsymbol{f}^{ext}_B\right) + \sum\left(\frac{\partial}{\partial\boldsymbol{\psi}}\left(C^{AB}\boldsymbol{m}_{0,B}\right) \delta\boldsymbol{\psi} + C^AB\delta\boldsymbol{m}_B^{ext}\right)\]

which gives the following stiffenning terms in the rigid-flex partition of the stiffness matrix:

\[K_{ss}^{M,R} \leftarrow +\sum\widetilde{C^{AB}_0 \boldsymbol{f}_{0,B}^{ext}}\]
\[K_{ss}^{M,\Psi} \leftarrow -\sum\widetilde{\boldsymbol{R}}\frac{\partial}{\partial\boldsymbol{\psi}} \left(C^{AB}\boldsymbol{f}_{0,B}\right)\]

and

\[K_{ss}^{M,\Psi} \leftarrow -\sum\frac{\partial}{\partial\boldsymbol{\psi}} \left(C^{AB}\boldsymbol{m}_{0,B}\right).\]
Parameters:

tsstr (sharpy.utils.datastructures.StructTimeStepInfo) – Linearisation time step.

linearise_gravity_forces(tsstr=None)[source]

Linearises gravity forces and includes the resulting terms in the C and K matrices. The method takes the linearisation condition (optional argument), linearises and updates:

  • Stiffness matrix

  • Damping matrix

  • Modal damping matrix

The method works for both the quaternion and euler angle orientation parametrisation.

Parameters:

tsstr (sharpy.utils.datastructures.StructTimeStepInfo) – Structural timestep at the linearisation point

Notes

The gravity forces are linearised to express them in terms of the beam formulation input variables:

  • Nodal forces: \(\delta \mathbf{f}_A\)

  • Nodal moments: \(\delta(T^T \mathbf{m}_B)\)

  • Total forces (rigid body equations): \(\delta \mathbf{F}_A\)

  • Total moments (rigid body equations): \(\delta \mathbf{M}_A\)

Gravity forces are naturally expressed in G (inertial) frame

\[\mathbf{f}_{G,0} = \mathbf{M\,g}\]

where the \(\mathbf{M}\) is the tangent mass matrix obtained at the linearisation reference.

To obtain the gravity forces expressed in A frame we make use of the projection matrices

\[\mathbf{f}_A = C^{AG}(\boldsymbol{\chi}) \mathbf{f}_{G,0}\]

that projects a vector in the inertial frame G onto the body attached frame A.

The projection of a vector can then be linearised as

\[\delta \mathbf{f}_A = C^{AG} \delta \mathbf{f}_{G,0} + \frac{\partial}{\partial \boldsymbol{\chi}}(C^{AG} \mathbf{f}_{G,0}) \delta\boldsymbol{\chi}.\]
  • Nodal forces:

    The linearisation of the gravity forces acting at each node is simply

    \[\delta \mathbf{f}_A = + \frac{\partial}{\partial \boldsymbol{\chi}}(C^{AG} \mathbf{f}_{G,0}) \delta\boldsymbol{\chi}\]

    where it is assumed that \(\delta\mathbf{f}_G = 0\).

  • Nodal moments:

    The gravity moments can be expressed in the local node frame of reference B by

    \[\mathbf{m}_B = \tilde{X}_{B,CG}C^{BA}(\Psi)C^{AG}(\boldsymbol{\chi})\mathbf{f}_{G,0}\]

    The linearisation is given by:

    \[\delta \mathbf{m}_B = \tilde{X}_{B,CG} \left(\frac{\partial}{\partial\Psi}(C^{BA}\mathbf{f}_{A,0})\delta\Psi + C^{BA}\frac{\partial}{\partial\boldsymbol{\chi}}(C^{AG}\mathbf{f}_{G,0})\delta\boldsymbol{\chi}\right)\]

    However, recall that the input moments are defined in tangential space \(\delta(T^\top\mathbf{m}_B)\) whose linearised expression is

    \[\delta(T^T(\Psi) \mathbf{m}_B) = T_0^T \delta \mathbf{m}_B + \frac{\partial}{\partial \Psi}(T^T \mathbf{m}_{B,0})\delta\Psi\]

    where the \(\delta \mathbf{m}_B\) term has been defined above.

  • Total forces:

    The total forces include the contribution from all flexible degrees of freedom as well as the gravity forces arising from the mass at the clamped node

    \[\mathbf{F}_A = \sum_n \mathbf{f}_A + \mathbf{f}_{A,clamped}\]

    which becomes

    \[\delta \mathbf{F}_A = \sum_n \delta \mathbf{f}_A + \frac{\partial}{\partial\boldsymbol{\chi}}\left(C^{AG}\mathbf{f}_{G,clamped}\right) \delta\boldsymbol{\chi}.\]
  • Total moments:

    The total moments, as opposed to the nodal moments, are expressed in A frame and again require the addition of the moments from the flexible structural nodes as well as the ones from the clamped node itself.

    \[\mathbf{M}_A = \sum_n \tilde{X}_{A,n}^{CG} C^{AG} \mathbf{f}_{n,G} + \tilde{X}_{A,clamped}C^{AG}\mathbf{f}_{G, clamped}\]

    where \(X_{A,n}^{CG} = R_{A,n} + C^{AB}(\Psi)X_{B,n}^{CG}\). Its linearised form is

    \[\delta X_{A,n}^{CG} = \delta R_{A,n} + \frac{\partial}{\partial \Psi}(C^{AB} X_{B,CG})\delta\Psi\]

    Therefore, the overall linearisation of the total moment is defined as

    \[\delta \mathbf{M}_A = \tilde{X}_{A,total}^{CG} \frac{\partial}{\partial \boldsymbol{\chi}}(C^{AG}\mathbf{F}_{G, total}) \delta \boldsymbol{\chi} -\sum_n \tilde{C}^{AG}\mathbf{f}_{G,0} \delta X_{A,n}^{CG}\]

    where \(X_{A, total}\) is the centre of gravity of the entire system expressed in A frame and \(\mathbf{F}_{G, total}\) are the gravity forces of the overall system in G frame, including the contributions from the clamped node.

The linearisation introduces damping and stiffening terms since the \(\delta\boldsymbol{\chi}\) and \(\delta\boldsymbol{\Psi}\) terms are found in the damping and stiffness matrices respectively.

Therefore, the beam matrices need updating to account for these terms:

  • Terms from the linearisation of the nodal moments will be assembled in the rows corresponding to moment equations and columns corresponding to the cartesian rotation vector

    \[K_{ss}^{m,\Psi} \leftarrow -T_0^T \tilde{X}_{B,CG} \frac{\partial}{\partial\Psi}(C^{BA}\mathbf{f}_{A,0}) -\frac{\partial}{\partial \Psi}(T^T \mathbf{m}_{B,0})\]
  • Terms from the linearisation of the translation forces with respect to the orientation are assembled in the damping matrix, the rows corresponding to translational forces and columns to orientation degrees of freedom

    \[C_{sr}^{f,\boldsymbol{\chi}} \leftarrow - \frac{\partial}{\partial \boldsymbol{\chi}}(C^{AG} \mathbf{f}_{G,0})\]
  • Terms from the linearisation of the moments with respect to the orientation are assembled in the damping matrix, with the rows correspondant to the moments and the columns to the orientation degrees of freedom

    \[C_{sr}^{m,\boldsymbol{\chi}} \leftarrow - T_0^T\tilde{X}_{B,CG}C^{BA}\frac{\partial}{\partial\boldsymbol{\chi}}(C^{AG}\mathbf{f}_{G,0})\]
  • Terms from the linearisation of the total forces with respect to the orientation correspond to the rigid body equations in the damping matrix, the rows to the translational forces and columns to the orientation

    \[C_{rr}^{F,\boldsymbol{\chi}} \leftarrow - \sum_n \frac{\partial}{\partial \boldsymbol{\chi}}(C^{AG} \mathbf{f}_{G,0})\]
  • Terms from the linearisation of the total moments with respect to the orientation correspond to the rigid body equations in the damping matrix, the rows to the moments and the columns to the orientation

    \[C_{rr}^{M,\boldsymbol{\chi}} \leftarrow - \sum_n\tilde{X}_{A,n}^{CG} \frac{\partial}{\partial \boldsymbol{\chi}}(C^{AG}\mathbf{f}_{G,0})\]
  • Terms from the linearisation of the total moments with respect to the nodal position \(R_A\) are included in the stiffness matrix, the rows corresponding to the moments in the rigid body equations and the columns to the nodal position

    \[K_{rs}^{M,R} \leftarrow + \sum_n \tilde{\mathbf{f}_{A,0}}\]
  • Terms from the linearisation of the total moments with respect to the cartesian rotation vector are included in the stiffness matrix, the rows corresponding to the moments in the rigid body equations and the columns to the cartesian rotation vector

    \[K_{rs}^{M, \Psi} \leftarrow + \sum_n \tilde{\mathbf{f}_{A,0}}\frac{\partial}{\partial \Psi}(C^{AB} X_{B,CG})\]
reshape_struct_input()[source]

Reshape structural input in a column vector

scale_system_normalised_time(time_ref)[source]

Scale the system with a normalised time step. The resulting time step is \(\Delta t = \Delta \bar{t}/t_{ref}\), where the over bar denotes dimensional time. The structural equations of motion are rescaled as:

\[\mathbf{M}\ddot{\boldsymbol{\eta}} + \mathbf{C} t_{ref} \dot{\boldsymbol{\eta}} + \mathbf{K} t_{ref}^2 \boldsymbol{\eta} = t_{ref}^2 \mathbf{N}\]

For aeroelastic applications, the reference time is usually defined using the semi-chord, \(b\), and the free stream velocity, \(U_\infty\).

\[t_{ref,ae} = \frac{b}{U_\infty}\]
Parameters:

time_ref (float) – Normalisation factor such that \(t/\bar{t}\) is non-dimensional.

tune_newmark_damp(amplification_factor=0.999)[source]

Tune artifical damping to achieve a percent reduction of the lower frequency (lower damped) mode

update_modal()[source]

Re-projects the full-states continuous-time structural dynamics equations

\[\mathbf{M}\,\mathbf{\ddot{x}} +\mathbf{C}\,\mathbf{\dot{x}} + \mathbf{K\,x} = \mathbf{F}\]

onto modal space. The modes used to project are controlled through the self.proj_modes={damped or undamped} attribute.

Warning

This method overrides SHARPy timestep_info results and requires Mstr, Cstr, Kstr to be available.

update_truncated_modes(nmodes)[source]

Updates the system to the specified number of modes

Parameters:

nmodes

Returns: