StepLinearUVLM
- class sharpy.solvers.steplinearuvlm.StepLinearUVLM[source]
Time domain aerodynamic solver that uses a linear UVLM formulation to be used with the
solvers.DynamicCoupled
solver.To use this solver, the
solver_id = StepLinearUVLM
must be given as the name for theaero_solver
is the case of an aeroelastic solver, where the setting below would be parsed throughaero_solver_settings
.Notes
The
integr_order
variable 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 whenintegr_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}\]If
track_body
isTrue
, the UVLM is projected onto a frameU
that is:Coincident with
G
at the linearisation timestep.Thence, rotates by the same quantity as the FoR
A
.
It is similar to a stability axes and is recommended any time rigid body dynamics are included.
See also
sharpy.sharpy.linear.assembler.linearuvlm.LinearUVLM
References
[1] 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
dt
float
Time step
0.1
integr_order
int
Integration order of the circulation derivative. Either
1
or2
.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
track_body
bool
UVLM inputs and outputs projected to coincide with lattice at linearisation
True
track_body_number
int
Frame of reference number to follow. If
-1
trackA
frame.-1
velocity_field_generator
str
Name of the velocity field generator to be used in the simulation
SteadyVelocityField
velocity_field_input
dict
Dictionary of settings for the velocity field generator
{}
vortex_radius
float
Distance between points below which induction is not computed
sharpy.utils.constants.vortex_radius_def
vortex_radius_wake_ind
float
Distance between points below which induction is not computed in the wake convection
sharpy.utils.constants.vortex_radius_def
cfl1
bool
If it is
True
, it assumes that the discretisation complies with CFL=1True
The settings that
ScalingDict
accepts are the following:Name
Type
Description
Default
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
- initialise(data, custom_settings=None, restart=False)[source]
Initialises the Linear UVLM aerodynamic solver and the chosen velocity generator.
Settings are parsed into the standard SHARPy settings format for solvers. It then checks whether there is any previous information about the linearised system (in order for a solution to be restarted without overwriting the linearisation).
If a linearised system does not exist, a linear UVLM system is created linearising about the current time step.
The reference values for the input and output are transformed into column vectors \(\mathbf{u}\) and \(\mathbf{y}\), respectively.
The information pertaining to the linear system is stored in a dictionary
self.data.aero.linear
within the maindata
variable.- Parameters:
data (PreSharpy) – class containing the problem information
custom_settings (dict) – custom settings dictionary
- pack_input_vector()[source]
Transform a SHARPy AeroTimestep instance into a column vector containing the input to the linear UVLM system.
\[[\zeta,\, \dot{\zeta}, u_{ext}] \longrightarrow \mathbf{u}\]If the
track_body
option is on, the function projects all the input into a frame that:is equal to the FoR G at time 0 (linearisation point)
rotates as the body frame specified in the
track_body_number
- Returns:
Input vector
- Return type:
np.ndarray
- static pack_state_vector(aero_tstep, aero_tstep_m1, dt, integr_order)[source]
Transform SHARPy Aerotimestep format into column vector containing the state information.
The state vector is of a different form depending on the order of integration chosen. If a second order scheme is chosen, the state includes the bound circulation at the previous timestep, hence the timestep information for the previous timestep shall be parsed.
The transformation is of the form:
If
integr_order==1
:\[\mathbf{x}_n = [\mathbf{\Gamma}^T_n,\, \mathbf{\Gamma_w}_n^T,\, \Delta t \,\mathbf{\dot{\Gamma}}_n^T]^T\]Else, if
integr_order==2
:\[\mathbf{x}_n = [\mathbf{\Gamma}_n^T,\, \mathbf{\Gamma_w}_n^T,\, \Delta t \,\mathbf{\dot{\Gamma}}_n^T,\, \mathbf{\Gamma}_{n-1}^T]^T\]
For the second order integration scheme, if the previous timestep information is not parsed, a first order stencil is employed to estimate the bound circulation at the previous timestep:
\[\mathbf{\Gamma}^{n-1} = \mathbf{\Gamma}^n - \Delta t \mathbf{\dot{\Gamma}}^n\]- Parameters:
aero_tstep (AeroTimeStepInfo) – Aerodynamic timestep information at the current timestep
n
.aero_tstep_m1 (AeroTimeStepInfo) –
- Returns:
State vector
- Return type:
np.ndarray
- run(**kwargs)[source]
Solve the linear aerodynamic UVLM model at the current time step
n
. The step increment is solved as:\[\begin{split}\mathbf{x}^n &= \mathbf{A\,x}^{n-1} + \mathbf{B\,u}^n \\ \mathbf{y}^n &= \mathbf{C\,x}^n + \mathbf{D\,u}^n\end{split}\]A change of state is possible in order to solve the system without the predictor term. In which case the system is solved by:
\[\begin{split}\mathbf{h}^n &= \mathbf{A\,h}^{n-1} + \mathbf{B\,u}^{n-1} \\ \mathbf{y}^n &= \mathbf{C\,h}^n + \mathbf{D\,u}^n\end{split}\]Variations are taken with respect to initial reference state. The state and input vectors for the linear UVLM system are of the form:
- If
integr_order==1
: - \[\mathbf{x}_n = [\delta\mathbf{\Gamma}^T_n,\, \delta\mathbf{\Gamma_w}_n^T,\, \Delta t \,\delta\mathbf{\dot{\Gamma}}_n^T]^T\]
- Else, if
integr_order==2
: - \[\mathbf{x}_n = [\delta\mathbf{\Gamma}_n^T,\, \delta\mathbf{\Gamma_w}_n^T,\, \Delta t \,\delta\mathbf{\dot{\Gamma}}_n^T,\, \delta\mathbf{\Gamma}_{n-1}^T]^T\]
- And the input vector:
- \[\mathbf{u}_n = [\delta\mathbf{\zeta}_n^T,\, \delta\dot{\mathbf{\zeta}}_n^T,\,\delta\mathbf{u_{ext}}^T_n]^T\]
where the subscript
n
refers to the time step.The linear UVLM system is then solved as detailed in
sharpy.linear.src.linuvlm.Dynamic.solve_step()
. The output is a column vector containing the aerodynamic forces at the panel vertices.To Do: option for impulsive start?
- Parameters:
aero_tstep (AeroTimeStepInfo) – object containing the aerodynamic data at the current time step
structure_tstep (StructTimeStepInfo) – object containing the structural data at the current time step
convect_wake (bool) – for backward compatibility only. The linear UVLM assumes a frozen wake geometry
dt (float) – time increment
t (float) – current time
unsteady_contribution (bool) – (backward compatibily). Unsteady aerodynamic effects are always included
- Returns:
updated
self.data
class with the new forces and circulation terms of the system- Return type:
PreSharpy
- If
- unpack_ss_vectors(y_n, x_n, u_n, aero_tstep)[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:equal to the FoR G at time 0 (linearisation point)
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 inx
,y
andz
. 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