Frequency
- class sharpy.linear.src.linuvlm.Frequency(tsdata, dt, integr_order=2, RemovePredictor=True, ScalingDict=None, UseSparse=True)[source]
Class for frequency description of linearised UVLM solution. Linearisation around steady-state are only supported. The class is built upon Static, and inherits all the methods contained there.
The class supports most of the features of Dynamics but has lower memory requirements of Dynamic, and should be preferred for:
producing memory and computationally cheap frequency responses
building reduced order models using RFA/polynomial fitting
Usage: Upon initialisation, the assemble method produces all the matrices required for the frequency description of the UVLM (see assemble for details). A state-space model is not allocated but:
Time stepping is also possible (but not implemented yet) as all the fundamental terms describing the UVLM equations are still produced (except the propagation of wake circulation)
ad-hoc methods for scaling, unscaling and frequency response are provided.
Input:
tsdata: aero timestep data from SHARPy solution
dt: time-step
integr_order=0,1,2: integration order for UVLM unsteady aerodynamic force. If 0, the derivative is computed exactly.
RemovePredictor=True: This flag is only used for the frequency response calculation. The frequency description, in fact, naturally arises without the predictor, but lags can be included during the frequency response calculation. See Dynamic documentation for more details.
ScalingDict=None: disctionary containing fundamental reference units
{'length': reference_length, 'speed': reference_speed, 'density': reference density}
used to derive scaling quantities for the state-space model variables. The scaling factors are stores in
self.ScalingFact.
Note that while time, circulation, angular speeds) are scaled accordingly, FORCES ARE NOT. These scale by qinf*b**2, where b is the reference length and qinf is the dinamic pressure.
UseSparse=False: builds the A and B matrices in sparse form. C and D are dense, hence the sparce format is not used.
- - nondimss
normalises matrices produced by the assemble method based on the scaling factors in self.ScalingFact.
- - dimss
inverse of nondimss.
- - assemble
builds matrices for UVLM minimal size description.
- - assemble_profiling
generate profiling report of the assembly and saves it into self.prof_out. To read the report:
import pstats p=pstats.Stats(self.prof_out)
- - freqresp
fast algorithm for frequency response.
Methods to implement:
solve_steady: runs freqresp at 0 frequency.
solve_step: solves one time-step
- assemble()[source]
Assembles matrices for minumal size frequency description of UVLM. The state equation is represented in the form:
\[\mathbf{A_0} \mathbf{\Gamma} + \mathbf{A_{w_0}} \mathbf{\Gamma_w} = \mathbf{B_0} \mathbf{u}\]While the output equation is as per the Dynamic class, namely:
\[\mathbf{y} = \mathbf{C} \mathbf{x} + \mathbf{D} \mathbf{u}\]where
\[\mathbf{x} = [\mathbf{\Gamma}; \mathbf{\Gamma_w}; \Delta\mathbf(\Gamma)]\]The propagation of wake circulation matrices are not produced as these are not required for frequency response analysis.
- assemble_profiling()[source]
Generate profiling report for assembly and save it in self.prof_out.
- To read the report:
import pstats p=pstats.Stats(self.prof_out)
- freqresp(kv, wake_prop_settings=None)[source]
Ad-hoc method for fast UVLM frequency response over the frequencies kv. The method, only requires inversion of a K x K matrix at each frequency as the equation for propagation of wake circulation are solved exactly.