LinearUVLM

class sharpy.linear.assembler.linearuvlm.LinearUVLM[source]

Linear UVLM System Assembler

Produces state-space model of the form

\[\begin{split}\mathbf{x}_{n+1} &= \mathbf{A}\,\mathbf{x}_n + \mathbf{B} \mathbf{u}_{n+1} \\ \mathbf{y}_n &= \mathbf{C}\,\mathbf{x}_n + \mathbf{D} \mathbf{u}_n\end{split}\]

where the state, inputs and outputs are:

\[\mathbf{x}_n = \{ \delta \mathbf{\Gamma}_n,\, \delta \mathbf{\Gamma_{w_n}},\, \Delta t\,\delta\mathbf{\Gamma}'_n,\, \delta\mathbf{\Gamma}_{n-1} \}\]
\[\mathbf{u}_n = \{ \delta\mathbf{\zeta}_n,\, \delta\mathbf{\zeta}'_n,\, \delta\mathbf{u}_{ext,n} \}\]
\[\mathbf{y} = \{\delta\mathbf{f}\}\]

with \(\mathbf{\Gamma}\in\mathbb{R}^{MN}\) being the vector of vortex circulations, \(\mathbf{\zeta}\in\mathbb{R}^{3(M+1)(N+1)}\) the vector of vortex lattice coordinates and \(\mathbf{f}\in\mathbb{R}^{3(M+1)(N+1)}\) the vector of aerodynamic forces and moments. Note that \((\bullet)'\) denotes a derivative with respect to time.

Note that the input is atypically defined at time n+1. If the setting remove_predictor = True the predictor term u_{n+1} is eliminated through the change of state[1]:

\[\begin{split}\mathbf{h}_n &= \mathbf{x}_n - \mathbf{B}\,\mathbf{u}_n \\\end{split}\]

such that:

\[\begin{split}\mathbf{h}_{n+1} &= \mathbf{A}\,\mathbf{h}_n + \mathbf{A\,B}\,\mathbf{u}_n \\ \mathbf{y}_n &= \mathbf{C\,h}_n + (\mathbf{C\,B}+\mathbf{D})\,\mathbf{u}_n\end{split}\]

which only modifies the equivalent \(\mathbf{B}\) and \(\mathbf{D}\) matrices.

The integr_order setting refers to the finite differencing scheme used to calculate the bound circulation derivative with respect to time \(\dot{\mathbf{\Gamma}}\). A first order scheme is used when integr_order == 1

\[\dot{\mathbf{\Gamma}}^{n+1} = \frac{\mathbf{\Gamma}^{n+1}-\mathbf{\Gamma}^n}{\Delta t}\]

If integr_order == 2 a higher order scheme is used (but it isn’t exactly second order accurate [1]).

\[\dot{\mathbf{\Gamma}}^{n+1} = \frac{3\mathbf{\Gamma}^{n+1}-4\mathbf{\Gamma}^n + \mathbf{\Gamma}^{n-1}} {2\Delta t}\]

References

[1] Franklin, GF and Powell, JD. Digital Control of Dynamic Systems, Addison-Wesley Publishing Company, 1980

[2] Maraniello, S., & Palacios, R.. State-Space Realizations and Internal Balancing in Potential-Flow Aerodynamics with Arbitrary Kinematics. AIAA Journal, 57(6), 1–14. 2019. https://doi.org/10.2514/1.J058153

The settings that this solver accepts are given by a dictionary, with the following key-value pairs:

Name Type Description Default Options
dt float Time step 0.1  
integr_order int Integration order of the circulation derivative. 2 1, 2
ScalingDict dict Dictionary of scaling factors to achieve normalised UVLM realisation. {}  
remove_predictor bool Remove the predictor term from the UVLM equations True  
use_sparse bool Assemble UVLM plant matrix in sparse format True  
density float Air density 1.225  
remove_inputs list(str) List of inputs to remove. u_gust to remove external velocity input. [] u_gust
gust_assembler str Selected linear gust assembler.   leading_edge
rom_method list(str) List of model reduction methods to reduce UVLM. []  
rom_method_settings dict Dictionary with settings for the desired ROM methods, where the name of the ROM method is the key to the dictionary {}  

The settings that this solver accepts are given by a dictionary, with the following key-value pairs:

Name Type Description Default Options
length float Reference length to be used for UVLM scaling 1.0  
speed float Reference speed to be used for UVLM scaling 1.0  
density float Reference density to be used for UVLM scaling 1.0  
assemble(track_body=False)[source]

Assembles the linearised UVLM system, removes the desired inputs and adds linearised control surfaces (if present).

With all possible inputs present, these are ordered as

\[\mathbf{u} = [\boldsymbol{\zeta},\,\dot{\boldsymbol{\zeta}},\,\mathbf{w},\,\delta]\]

Control surface inputs are ordered last as:

\[[\delta_1, \delta_2, \dots, \dot{\delta}_1, \dot{\delta_2}]\]
remove_inputs(remove_list=<class 'list'>)[source]

Remove certain inputs from the input vector

To do:
  • Support for block UVLM
Parameters:remove_list (list) – Inputs to remove
unpack_input_vector(u_n)[source]

Unpacks the input vector into the corresponding grid coordinates, velocities and external velocities.

Parameters:u_n (np.ndarray) – UVLM input vector. May contain control surface deflections and external velocities.
Returns:Tuple containing zeta, zeta_dot and u_ext, accounting for the effect of control surfaces.
Return type:tuple
unpack_ss_vector(data, x_n, aero_tstep, track_body=False)[source]

Transform column vectors used in the state space formulation into SHARPy format

The column vectors are transformed into lists with one entry per aerodynamic surface. Each entry contains a matrix with the quantities at each grid vertex.

\[\mathbf{y}_n \longrightarrow \mathbf{f}_{aero}\]
\[\mathbf{x}_n \longrightarrow \mathbf{\Gamma}_n,\, \mathbf{\Gamma_w}_n,\, \mathbf{\dot{\Gamma}}_n\]

If the track_body option is on, the output forces are projected from the linearization frame, to the G frame. Note that the linearisation frame is:

  1. equal to the FoR G at time 0 (linearisation point)
  2. rotates as the body frame specified in the track_body_number
Parameters:
  • y_n (np.ndarray) – Column output vector of linear UVLM system
  • x_n (np.ndarray) – Column state vector of linear UVLM system
  • u_n (np.ndarray) – Column input vector of linear UVLM system
  • aero_tstep (AeroTimeStepInfo) – aerodynamic timestep information class instance
Returns:

Tuple containing:

forces (list):

Aerodynamic forces in a list with n_surf entries. Each entry is a (6, M+1, N+1) matrix, where the first 3 indices correspond to the components in x, y and z. The latter 3 are zero.

gamma (list):

Bound circulation list with n_surf entries. Circulation is stored in an (M+1, N+1) matrix, corresponding to the panel vertices.

gamma_dot (list):

Bound circulation derivative list with n_surf entries. Circulation derivative is stored in an (M+1, N+1) matrix, corresponding to the panel vertices.

gamma_star (list):

Wake (free) circulation list with n_surf entries. Wake circulation is stored in an (M_star+1, N+1) matrix, corresponding to the panel vertices of the wake.

Return type:

tuple