def duhamel_integral_kussner(reduced_time: np.ndarray, gust_velocity: np.ndarray, velocity: float): """ Calculates the duhamel superposition integral of Kussner's problem. Given some arbitrary transverse velocity profile, the lift coefficient of a flat plate can be computed using this function Args: reduced_time (float,np.ndarray) : Reduced time, equal to the number of semichords travelled. See function reduced_time gust_velocity (np.ndarray) : The transverse velocity profile that the flate plate experiences velocity (float) :The velocity by which the flat plate enters the gust Returns: lift_coefficient (np.ndarray) : The lift coefficient history of the flat plate """ assert np.size(reduced_time) == np.size( gust_velocity ), "The velocity history and time must have the same length" dw_ds = np.gradient(gust_velocity) lift_coefficient = np.zeros_like(reduced_time) kussner = kussners_function(reduced_time) ds = np.gradient(reduced_time) for i, s in enumerate(reduced_time): lift_coefficient[i] = 2 * np.pi / velocity * ( gust_velocity[0] * kussner[i] + np.sum([dw_ds[j] * kussner[i - j] * ds[j] for j in range(i)])) return lift_coefficient
def calculate_reduced_time(time: Union[float, np.ndarray], velocity: Union[float, np.ndarray], chord: float) -> Union[float, np.ndarray]: """ Calculates reduced time from time in seconds and velocity history in m/s. For constant velocity it reduces to s = 2*U*t/c The reduced time is the number of semichords travelled by the airfoil/aircaft i.e. 2 / chord * integral from t0 to t of velocity dt Args: time (float,np.ndarray) : Time in seconds velocity (float,np.ndarray): Either a constant velocity or array of velocities at corresponding reduced times chord (float) : The chord of the airfoil Returns: The reduced time as an ndarray or float similar to the input. The first element is 0. """ if type(velocity) == float or type(velocity) == int: return 2 * velocity * time / chord else: assert np.size(velocity) == np.size( time), "The velocity history and time must have the same length" reduced_time = np.zeros_like(time) for i in range(len(time) - 1): reduced_time[i + 1] = reduced_time[i] + ( velocity[i + 1] + velocity[i]) / 2 * (time[i + 1] - time[i]) return 2 / chord * reduced_time
def duhamel_integral_wagner(reduced_time: np.ndarray, angle_of_attack: np.ndarray): """ Calculates the duhamel superposition integral of Wagner's problem. Given some arbitrary pitching profile, the lift coefficient of a flat plate can be computed using this function Args: reduced_time (float,np.ndarray) : Reduced time, equal to the number of semichords travelled. See function reduced_time angle_of_attack (np.ndarray) : The angle of attack as a function of reduced time of the flat plate Returns: lift_coefficient (np.ndarray) : The lift coefficient history of the flat plate """ assert np.size(reduced_time) == np.size( angle_of_attack ), "The pitching history and time must have the same length" da_ds = np.gradient(angle_of_attack) lift_coefficient = np.zeros_like(reduced_time) wagner = wagners_function(reduced_time) ds = np.gradient(reduced_time) for i, s in enumerate(reduced_time): lift_coefficient[i] = 2 * np.pi * ( angle_of_attack[0] * wagner[i] + np.sum([da_ds[j] * wagner[i - j] * ds[j] for j in range(i)])) return lift_coefficient
def calculate_lift_due_to_pitching_profile( reduced_time: np.ndarray, angle_of_attack: Union[Callable[[float], float], float] # In degrees ): """ Calculates the duhamel superposition integral of Wagner's problem. Given some arbitrary pitching profile. The lift coefficient as a function of reduced time of a flat plate can be computed using this function Args: reduced_time (float,np.ndarray) : Reduced time, equal to the number of semichords travelled. See function reduced_time angle_of_attack (Callable[[float],float]) : The angle of attack as a function of reduced time of the flat plate. Must be a Callable that takes reduced time and returns angle of attack Returns: lift_coefficient (np.ndarray) : The lift coefficient history of the flat plate """ assert (reduced_time >= 0).all(), "Please use positive time. Negative time not supported" if isinstance(angle_of_attack, float) or isinstance(angle_of_attack, int): def AoA_function(reduced_time): return np.deg2rad(angle_of_attack) else: def AoA_function(reduced_time): return np.deg2rad(angle_of_attack(reduced_time)) def dW_ds(reduced_time): return (0.1005 * np.exp(-0.3 * reduced_time) + 0.00750075 * np.exp(-0.0455 * reduced_time)) def integrand(sigma, s): if dW_ds(sigma) < 0: dW_ds(sigma) return dW_ds(sigma) * AoA_function(s - sigma) lift_coefficient = np.zeros_like(reduced_time) for i, s in enumerate(reduced_time): I = quad(integrand, 0, s, args=s)[0] # print(I) lift_coefficient[i] = 2 * np.pi * ( AoA_function(s) * wagners_function(0) + I) return lift_coefficient
def calculate_lift_due_to_transverse_gust( reduced_time: np.ndarray, gust_velocity_profile: Callable[[float], float], plate_velocity: float, angle_of_attack: Union[float, Callable[[float], float]] = 0, # In Degrees chord: float = 1): """ Calculates the lift (as a function of reduced time) caused by an arbitrary transverse gust profile by computing duhamel superposition integral of Kussner's problem at a constant angle of attack Args: reduced_time (float,np.ndarray) : Reduced time, equal to the number of semichords travelled. See function reduced_time gust_velocity_profile (Callable[[float],float]) : The transverse velocity profile that the flate plate experiences. Must be a function that takes reduced time and returns a velocity plate_velocity (float) :The velocity by which the flat plate enters the gust angle_of_attack (Union[float,Callable[[float],float]]) : The angle of attack, in degrees. Can either be a float for constant angle of attack or a Callable that takes reduced time and returns angle of attack chord (float) : The chord of the plate in meters Returns: lift_coefficient (np.ndarray) : The lift coefficient history of the flat plate """ assert type( angle_of_attack ) != np.ndarray, "Please provide either a Callable or a float for the angle of attack" if isinstance(angle_of_attack, float) or isinstance(angle_of_attack, int): def AoA_function(reduced_time): return np.deg2rad(angle_of_attack) else: def AoA_function(reduced_time): return np.deg2rad(angle_of_attack(reduced_time)) def dK_ds(reduced_time): return (0.065 * np.exp(-0.13 * reduced_time) + 0.5 * np.exp(-reduced_time)) def integrand(sigma, s, chord): offset = chord / 2 * (1 - np.cos(AoA_function(s - sigma))) return (dK_ds(sigma) * gust_velocity_profile(s - sigma - offset) * np.cos(AoA_function(s - sigma))) lift_coefficient = np.zeros_like(reduced_time) for i, s in enumerate(reduced_time): I = quad(integrand, 0, s, args=(s, chord))[0] lift_coefficient[i] = 2 * np.pi * I / plate_velocity return lift_coefficient
def airfoil_coefficients_post_stall( airfoil: Airfoil, alpha: float, ): """ Estimates post-stall aerodynamics of an airfoil. Uses methods given in: Truong, V. K. "An analytical model for airfoil aerodynamic characteristics over the entire 360deg angle of attack range". J. Renewable Sustainable Energy. 2020. doi: 10.1063/1.5126055 Args: airfoil: op_point: Returns: """ sina = np.sind(alpha) cosa = np.cosd(alpha) ##### Normal force calulation # Cd90_fp = aerolib.Cd_flat_plate_normal() # TODO implement # Cd90_0 = Cd90_fp - 0.83 * airfoil.LE_radius() - 1.46 / 2 * airfoil.max_thickness() + 1.46 * airfoil.max_camber() # Cd270_0 = Cd90_fp - 0.83 * airfoil.LE_radius() - 1.46 / 2 * airfoil.max_thickness() - 1.46 * airfoil.max_camber() ### Values for NACA0012 Cd90_0 = 2.08 pn2_star = 8.36e-2 pn3_star = 4.06e-1 pt1_star = 9.00e-2 pt2_star = -1.78e-1 pt3_star = -2.98e-1 Cd90 = Cd90_0 + pn2_star * cosa + pn3_star * cosa**2 CN = Cd90 * sina ##### Tangential force calculation CT = (pt1_star + pt2_star * cosa + pt3_star * cosa**3) * sina**2 ##### Conversion to wind axes CL = CN * cosa + CT * sina CD = CN * sina - CT * cosa CM = np.zeros_like(CL) # TODO return CL, CD, CM
def wing_aerodynamics( self, wing: Wing, ): """ Estimates the aerodynamic forces, moments, and derivatives on a wing in isolation. Moments are given with the reference at Wing [0, 0, 0]. Args: wing: A Wing object that you wish to analyze. op_point: The OperatingPoint that you wish to analyze the fuselage at. TODO account for wing airfoil pitching moment Returns: """ ##### Alias a few things for convenience op_point = self.op_point wing_options = self.get_options(wing) ##### Compute general wing properties and things to be used in sectional analysis. sweep = wing.mean_sweep_angle() AR = wing.aspect_ratio() mach = op_point.mach() q = op_point.dynamic_pressure() CL_over_Cl = aerolib.CL_over_Cl(aspect_ratio=AR, mach=mach, sweep=sweep) oswalds_efficiency = aerolib.oswalds_efficiency( taper_ratio=wing.taper_ratio(), aspect_ratio=AR, sweep=sweep, fuselage_diameter_to_span_ratio=0 # an assumption ) areas = wing.area(_sectional=True) aerodynamic_centers = wing.aerodynamic_center(_sectional=True) F_g = [0, 0, 0] M_g = [0, 0, 0] ##### Iterate through the wing sections. for sect_id in range(len(wing.xsecs) - 1): ##### Identify the wing cross sections adjacent to this wing section. xsec_a = wing.xsecs[sect_id] xsec_b = wing.xsecs[sect_id + 1] ##### When linearly interpolating, weight things by the relative chord. a_weight = xsec_a.chord / (xsec_a.chord + xsec_b.chord) b_weight = xsec_b.chord / (xsec_a.chord + xsec_b.chord) ##### Compute the local frame of this section, and put the z (normal) component into wind axes. xg_local, yg_local, zg_local = wing._compute_frame_of_section( sect_id) sect_aerodynamic_center = aerodynamic_centers[sect_id] sect_z_w = op_point.convert_axes( x_from=zg_local[0], y_from=zg_local[1], z_from=zg_local[2], from_axes="geometry", to_axes="wind", ) ##### Compute the generalized angle of attack, so the geometric alpha that the wing section "sees". velocity_vector_b_from_freestream = op_point.convert_axes( x_from=-op_point.velocity, y_from=0, z_from=0, from_axes="wind", to_axes="body") velocity_vector_b_from_rotation = np.cross(op_point.convert_axes( sect_aerodynamic_center[0], sect_aerodynamic_center[1], sect_aerodynamic_center[2], from_axes="geometry", to_axes="body"), [op_point.p, op_point.q, op_point.r], manual=True) velocity_vector_b = [ velocity_vector_b_from_freestream[i] + velocity_vector_b_from_rotation[i] for i in range(3) ] velocity_mag_b = np.sqrt( sum([comp**2 for comp in velocity_vector_b])) velocity_dir_b = [ velocity_vector_b[i] / velocity_mag_b for i in range(3) ] sect_z_b = op_point.convert_axes( x_from=zg_local[0], y_from=zg_local[1], z_from=zg_local[2], from_axes="geometry", to_axes="body", ) vel_dot_normal = np.dot(velocity_dir_b, sect_z_b, manual=True) sect_alpha_generalized = 90 - np.arccosd(vel_dot_normal) def get_deflection(xsec): n_surfs = len(xsec.control_surfaces) if n_surfs == 0: return 0 elif n_surfs == 1: surf = xsec.control_surfaces[0] return surf.deflection else: raise NotImplementedError( "AeroBuildup currently cannot handle multiple control surfaces attached to a given WingXSec." ) ##### Compute sectional lift at cross sections using lookup functions. Merge them linearly to get section CL. xsec_a_Cl_incompressible = xsec_a.airfoil.CL_function( alpha=sect_alpha_generalized, Re=op_point.reynolds(xsec_a.chord), mach= 0, # Note: this is correct, mach correction happens in 2D -> 3D step deflection=get_deflection(xsec_a)) xsec_b_Cl_incompressible = xsec_b.airfoil.CL_function( alpha=sect_alpha_generalized, Re=op_point.reynolds(xsec_b.chord), mach= 0, # Note: this is correct, mach correction happens in 2D -> 3D step deflection=get_deflection(xsec_b)) sect_CL = (xsec_a_Cl_incompressible * a_weight + xsec_b_Cl_incompressible * b_weight) * CL_over_Cl ##### Compute sectional drag at cross sections using lookup functions. Merge them linearly to get section CD. xsec_a_Cd_profile = xsec_a.airfoil.CD_function( alpha=sect_alpha_generalized, Re=op_point.reynolds(xsec_a.chord), mach=mach, deflection=get_deflection(xsec_a)) xsec_b_Cd_profile = xsec_b.airfoil.CD_function( alpha=sect_alpha_generalized, Re=op_point.reynolds(xsec_b.chord), mach=mach, deflection=get_deflection(xsec_b)) sect_CDp = (xsec_a_Cd_profile * a_weight + xsec_b_Cd_profile * b_weight) ##### Compute induced drag from local CL and full-wing properties (AR, e) sect_CDi = (sect_CL**2 / (np.pi * AR * oswalds_efficiency)) ##### Total the drag. sect_CD = sect_CDp + sect_CDi ##### Go to dimensional quantities using the area. area = areas[sect_id] sect_L = q * area * sect_CL sect_D = q * area * sect_CD ##### Compute the direction of the lift by projecting the section's normal vector into the plane orthogonal to the freestream. sect_L_direction_w = (np.zeros_like(sect_z_w[0]), sect_z_w[1] / np.sqrt(sect_z_w[1]**2 + sect_z_w[2]**2), sect_z_w[2] / np.sqrt(sect_z_w[1]**2 + sect_z_w[2]**2)) sect_L_direction_g = op_point.convert_axes(*sect_L_direction_w, from_axes="wind", to_axes="geometry") ##### Compute the direction of the drag by aligning the drag vector with the freestream vector. sect_D_direction_w = (-1, 0, 0) sect_D_direction_g = op_point.convert_axes(*sect_D_direction_w, from_axes="wind", to_axes="geometry") ##### Compute the force vector in geometry axes. sect_F_g = [ sect_L * sect_L_direction_g[i] + sect_D * sect_D_direction_g[i] for i in range(3) ] ##### Compute the moment vector in geometry axes. sect_M_g = np.cross(sect_aerodynamic_center, sect_F_g, manual=True) ##### Add section forces and moments to overall forces and moments F_g = [F_g[i] + sect_F_g[i] for i in range(3)] M_g = [M_g[i] + sect_M_g[i] for i in range(3)] ##### Treat symmetry if wing.symmetric: ##### Compute the local frame of this section, and put the z (normal) component into wind axes. sym_sect_aerodynamic_center = aerodynamic_centers[sect_id] sym_sect_aerodynamic_center[1] *= -1 sym_sect_z_w = op_point.convert_axes( x_from=zg_local[0], y_from=-zg_local[1], z_from=zg_local[2], from_axes="geometry", to_axes="wind", ) ##### Compute the generalized angle of attack, so the geometric alpha that the wing section "sees". sym_velocity_vector_b_from_freestream = op_point.convert_axes( x_from=-op_point.velocity, y_from=0, z_from=0, from_axes="wind", to_axes="body") sym_velocity_vector_b_from_rotation = np.cross( op_point.convert_axes(sym_sect_aerodynamic_center[0], sym_sect_aerodynamic_center[1], sym_sect_aerodynamic_center[2], from_axes="geometry", to_axes="body"), [op_point.p, op_point.q, op_point.r], manual=True) sym_velocity_vector_b = [ sym_velocity_vector_b_from_freestream[i] + sym_velocity_vector_b_from_rotation[i] for i in range(3) ] sym_velocity_mag_b = np.sqrt( sum([comp**2 for comp in sym_velocity_vector_b])) sym_velocity_dir_b = [ sym_velocity_vector_b[i] / sym_velocity_mag_b for i in range(3) ] sym_sect_z_b = op_point.convert_axes( x_from=zg_local[0], y_from=-zg_local[1], z_from=zg_local[2], from_axes="geometry", to_axes="body", ) sym_vel_dot_normal = np.dot(sym_velocity_dir_b, sym_sect_z_b, manual=True) sym_sect_alpha_generalized = 90 - np.arccosd( sym_vel_dot_normal) def get_deflection(xsec): n_surfs = len(xsec.control_surfaces) if n_surfs == 0: return 0 elif n_surfs == 1: surf = xsec.control_surfaces[0] return surf.deflection if surf.symmetric else -surf.deflection else: raise NotImplementedError( "AeroBuildup currently cannot handle multiple control surfaces attached to a given WingXSec." ) ##### Compute sectional lift at cross sections using lookup functions. Merge them linearly to get section CL. sym_xsec_a_Cl_incompressible = xsec_a.airfoil.CL_function( alpha=sym_sect_alpha_generalized, Re=op_point.reynolds(xsec_a.chord), mach= 0, # Note: this is correct, mach correction happens in 2D -> 3D step deflection=get_deflection(xsec_a)) sym_xsec_b_Cl_incompressible = xsec_b.airfoil.CL_function( alpha=sym_sect_alpha_generalized, Re=op_point.reynolds(xsec_b.chord), mach= 0, # Note: this is correct, mach correction happens in 2D -> 3D step deflection=get_deflection(xsec_b)) sym_sect_CL = ( sym_xsec_a_Cl_incompressible * a_weight + sym_xsec_b_Cl_incompressible * b_weight) * CL_over_Cl ##### Compute sectional drag at cross sections using lookup functions. Merge them linearly to get section CD. sym_xsec_a_Cd_profile = xsec_a.airfoil.CD_function( alpha=sym_sect_alpha_generalized, Re=op_point.reynolds(xsec_a.chord), mach=mach, deflection=get_deflection(xsec_a)) sym_xsec_b_Cd_profile = xsec_b.airfoil.CD_function( alpha=sym_sect_alpha_generalized, Re=op_point.reynolds(xsec_b.chord), mach=mach, deflection=get_deflection(xsec_b)) sym_sect_CDp = (sym_xsec_a_Cd_profile * a_weight + sym_xsec_b_Cd_profile * b_weight) ##### Compute induced drag from local CL and full-wing properties (AR, e) sym_sect_CDi = (sym_sect_CL**2 / (np.pi * AR * oswalds_efficiency)) ##### Total the drag. sym_sect_CD = sym_sect_CDp + sym_sect_CDi ##### Go to dimensional quantities using the area. area = areas[sect_id] sym_sect_L = q * area * sym_sect_CL sym_sect_D = q * area * sym_sect_CD ##### Compute the direction of the lift by projecting the section's normal vector into the plane orthogonal to the freestream. sym_sect_L_direction_w = ( np.zeros_like(sym_sect_z_w[0]), sym_sect_z_w[1] / np.sqrt(sym_sect_z_w[1]**2 + sym_sect_z_w[2]**2), sym_sect_z_w[2] / np.sqrt(sym_sect_z_w[1]**2 + sym_sect_z_w[2]**2)) sym_sect_L_direction_g = op_point.convert_axes( *sym_sect_L_direction_w, from_axes="wind", to_axes="geometry") ##### Compute the direction of the drag by aligning the drag vector with the freestream vector. sym_sect_D_direction_w = (-1, 0, 0) sym_sect_D_direction_g = op_point.convert_axes( *sym_sect_D_direction_w, from_axes="wind", to_axes="geometry") ##### Compute the force vector in geometry axes. sym_sect_F_g = [ sym_sect_L * sym_sect_L_direction_g[i] + sym_sect_D * sym_sect_D_direction_g[i] for i in range(3) ] ##### Compute the moment vector in geometry axes. sym_sect_M_g = np.cross(sym_sect_aerodynamic_center, sym_sect_F_g, manual=True) ##### Add section forces and moments to overall forces and moments F_g = [F_g[i] + sym_sect_F_g[i] for i in range(3)] M_g = [M_g[i] + sym_sect_M_g[i] for i in range(3)] ##### Convert F_g and M_g to body and wind axes for reporting. F_b = op_point.convert_axes(*F_g, from_axes="geometry", to_axes="body") F_w = op_point.convert_axes(*F_b, from_axes="body", to_axes="wind") M_b = op_point.convert_axes(*M_g, from_axes="geometry", to_axes="body") M_w = op_point.convert_axes(*M_b, from_axes="body", to_axes="wind") return { "F_g": F_g, "F_b": F_b, "F_w": F_w, "M_g": M_g, "M_b": M_b, "M_w": M_w, "L": -F_w[2], "Y": F_w[1], "D": -F_w[0], "l_b": M_b[0], "m_b": M_b[1], "n_b": M_b[2] }
def constrain_derivative( self, derivative: cas.MX, variable: cas.MX, with_respect_to: Union[np.ndarray, cas.MX], method: str = "midpoint", ) -> None: """ Adds a constraint to the optimization problem such that: d(variable) / d(with_respect_to) == derivative Can be used directly; also called indirectly by opti.derivative_of() for implicit derivative creation. Args: derivative: The derivative that is to be constrained here. variable: The variable or quantity that you are taking the derivative of. The "numerator" of the derivative, in colloquial parlance. with_respect_to: The variable or quantity that you are taking the derivative with respect to. The "denominator" of the derivative, in colloquial parlance. In a typical example case, this `with_respect_to` parameter would be time. Please make sure that the value of this parameter is monotonically increasing, otherwise you may get nonsensical answers. method: The type of integrator to use to define this derivative. Options are: * "forward euler" - a first-order-accurate forward Euler method Citation: https://en.wikipedia.org/wiki/Euler_method * "backwards euler" - a first-order-accurate backwards Euler method Citation: https://en.wikipedia.org/wiki/Backward_Euler_method * "midpoint" or "trapezoid" - a second-order-accurate midpoint method Citation: https://en.wikipedia.org/wiki/Midpoint_method * "simpson" - Simpson's rule for integration Citation: https://en.wikipedia.org/wiki/Simpson%27s_rule * "runge-kutta" or "rk4" - a fourth-order-accurate Runge-Kutta method. I suppose that technically, "forward euler", "backward euler", and "midpoint" are all (lower-order) Runge-Kutta methods... Citation: https://en.wikipedia.org/wiki/Runge%E2%80%93Kutta_methods#The_Runge%E2%80%93Kutta_method * "runge-kutta-3/8" - A modified version of the Runge-Kutta 4 proposed by Kutta in 1901. Also fourth-order-accurate, but all of the error coefficients are smaller than they are in the standard Runge-Kutta 4 method. The downside is that more floating point operations are required per timestep, as the Butcher tableau is more dense (i.e. not banded). Citation: Kutta, Martin (1901), "Beitrag zur näherungsweisen Integration totaler Differentialgleichungen", Zeitschrift für Mathematik und Physik, 46: 435–453 Note that all methods are expressed as integrators rather than differentiators; this prevents singularities from forming in the limit of timestep approaching zero. (For those coming from the PDE world, this is analogous to using finite volume methods rather than finite difference methods to allow shock capturing.) Returns: None (adds constraint in-place). """ try: d_var = np.diff(variable) except ValueError: d_var = np.diff(np.zeros_like(with_respect_to)) try: derivative[0] except (TypeError, IndexError): derivative = np.full_like(with_respect_to, fill_value=derivative) d_time = np.diff(with_respect_to) # Calculate the timestep # TODO scale constraints by variable scale? # TODO make if method == "forward euler" or method == "forward" or method == "forwards": # raise NotImplementedError self.subject_to(d_var == derivative[:-1] * d_time) elif method == "backward euler" or method == "backward" or method == "backwards": # raise NotImplementedError self.subject_to(d_var == derivative[1:] * d_time) elif method == "midpoint" or method == "trapezoid" or method == "trapezoidal": self.subject_to(d_var == np.trapz(derivative) * d_time, ) elif method == "simpson": raise NotImplementedError elif method == "runge-kutta" or method == "rk4": raise NotImplementedError elif method == "runge-kutta-3/8": raise NotImplementedError else: raise ValueError("Bad value of `method`!")
angle_of_attack[0] * wagner[i] + np.sum([da_ds[j] * wagner[i - j] * ds[j] for j in range(i)])) return lift_coefficient if __name__ == "__main__": n = 1000 n1 = int(n / 3) n2 = int(2 * n / 3) time = np.linspace(0, 100, n) velocity = 0.2 chord = 1 reduced_time = calculate_reduced_time(time, velocity, chord) gust_velocity = np.zeros_like(reduced_time) gust_velocity[n1:n2] = velocity angle_of_attack = 20 * np.deg2rad(np.sin(reduced_time)) #angle_of_attack[n1:n2] = np.deg2rad(-20) cl_k = duhamel_integral_kussner(reduced_time, gust_velocity, velocity) cl_w = duhamel_integral_wagner(reduced_time, angle_of_attack) plt.figure(dpi=300) plt.plot(reduced_time, cl_w, label="wagner") plt.plot(reduced_time, cl_k, label="kussner") plt.plot(reduced_time, cl_k + cl_w, label="total") plt.xlabel("Reduced time, t*") plt.ylabel("$C_\ell$") plt.legend()