AcousticAnalogies.AcousticObserver
— TypeSupertype for an object that recieves a noise prediction when combined with an acoustic analogy source; computational equivalent of a microphone.
(obs::AcousticObserver)(t)
Calculate the position of the acoustic observer at time t
.
AcousticAnalogies.CompactSourceElement
— MethodCompactSourceElement(ρ0, c0, r, θ, Δr, Λ, fn, fc, τ)
Construct a source element to be used with the compact form of Farassat's formulation 1A.
Arguments
- ρ0: Ambient air density (kg/m^3)
- c0: Ambient speed of sound (m/s)
- r: radial coordinate of the element in the blade-fixed coordinate system (m)
- θ: angular offest of the element in the blade-fixed coordinate system (rad)
- Δr: length of the element (m)
- Λ: cross-sectional area of the element (m^2)
- fn: normal load per unit span on the fluid (N/m)
- fr: radial load on the fluid (N/m)
- fc: circumferential load on the fluid (N/m)
- τ: source time (s)
AcousticAnalogies.CompactSourceElement
— MethodCompactSourceElement(rotor::CCBlade.Rotor, section::CCBlade.Section, op::CCBlade.OperatingPoint, out::CCBlade.Outputs, θ, Δr, area_per_chord2, τ)
Construct a source element to be used with the compact form of Farassat's formulation 1A from CCBlade objects.
The source element's position is calculated from section.r
, rotor.precone
, and the θ
argument using
sθ, cθ = sincos(θ)
spc, cpc = sincos(precone)
y0dot = [r*spc, r*cpc*cθ, r*cpc*sθ]
where y0dot
is the position of the source element.
Arguments
rotor::CCBlade.Rotor
: CCBlade rotor object, needed for the precone angle.precone.section::CCBlade.Section
: CCBlade section object, needed for the radial location and chord length of the element.op::CCBlade.OperatingPoint
: CCBlade operating point, needed for atmospheric properties.out::CCBlade.Outputs
: CCBlade outputs object, needed for the loading.θ
: polar coordinate of the element, in radians.Δr
: length of the element.area_per_chord2
: cross-sectional area divided by the chord squared of the element.τ
: source time of the element.
AcousticAnalogies.ConstVelocityAcousticObserver
— TypeConstVelocityAcousticObserver(t0, x0, v)
Construct an acoustic observer moving with a constant velocity v
, located at x0
at time t0
.
AcousticAnalogies.F1AOutput
— TypeOutput of the F1A calculation: the acoustic pressure value at time t
, broken into monopole component p_m
and dipole component p_d
.
AcousticAnalogies.F1APressureTimeHistory
— TypeF1APressureTimeHistory(apth::AbstractArray{<:F1AOutput}, period::AbstractFloat, n::Integer, axis::Integer=1)
Construct an F1APressureTimeHistory
struct
suitable for containing an acoustic prediction from an array of F1AOutput
struct
.
The elapsed time and length of the returned F1APressureTimeHistory
will be period
and n
, respectively. axis
indicates which axis the apth
struct
s time varies. (period
, n
, axis
are passed to common_obs_time
.)
AcousticAnalogies.F1APressureTimeHistory
— MethodF1APressureTimeHistory([T=Float64,] n, dt, t0)
Construct an F1APressureTimeHistory
struct
suitable for containing an acoustic prediction of length n
, starting at time t0
with time step dt
.
AcousticAnalogies.StationaryAcousticObserver
— TypeStationaryAcousticObserver(x)
Construct an acoustic observer that does not move with position x
(m).
KinematicCoordinateTransformations.KinematicTransformation
— Method(trans::KinematicTransformation)(se::CompactSourceElement)
Transform the position and forces of a source element according to the coordinate system transformation trans
.
AcousticAnalogies.adv_time
— Methodadv_time(se::CompactSourceElement, obs::AcousticObserver)
Calculate the time an acoustic wave emmited by source se
at time se.τ
is recieved by observer obs
.
AcousticAnalogies.combine
— Functioncombine(apth::AbstractArray{<:F1AOutput}, period::AbstractFloat, n::Integer, axis=1; f_interp=akima)
Combine the acoustic pressures of multiple sources (apth
) into a single acoustic pressure time history on a time grid of size n
extending over time length period
.
AcousticAnalogies.combine!
— Methodcombine!(apth_out::F1APressureTimeHistory, apth::AbstractArray{<:F1AOutput}, axis; f_interp=akima)
Combine the acoustic pressures of multiple sources (apth
) into a single acoustic pressure time history apth_out
.
The input acoustic pressures apth
are interpolated onto the time grid returned by time(apth_out)
. The interpolation is performed by the function f_intep(xpt, ypt, x)
, where xpt
and ytp
are the input grid and function values, respectively, and x
is the output grid.
AcousticAnalogies.common_obs_time
— Functioncommon_obs_time(apth::AbstractArray{<:F1AOutput}, period, n, axis=1)
Return a suitable time range for the collection of F1A acoustic pressures in apth
.
The time range will begin near the latest start time of the acoustic pressures in apth
, and be an AbstractVector
(really a StepRangeLen
) of size n
and of time length period
. axis
indicates which axis of apth
the time for a source varies.
AcousticAnalogies.endpoints
— Methodendpoints(se::CompactSourceElement)
Return the Tuple containing the endpoint locations of the compact source element se
.
AcousticAnalogies.f1a
— MethodAcousticAnalogies.f1a
— Methodf1a(se::CompactSourceElement, obs::AcousticObserver)
Calculate the acoustic pressure emitted by source element se
and recieved by observer obs
, returning an F1AOutput
object.
AcousticAnalogies.get_ccblade_dradii
— Methodget_ccblade_dradii(rotor::CCBlade.Rotor, sections::Vector{CCBlade.Section})
Construct and return a Vector of the lengths of each CCBlade section.
AcousticAnalogies.get_dradii
— Methodget_dradii(radii, Rhub, Rtip)
Compute the spacing between blade elements given the radial locations of the element midpoints in radii
and the hub and tip radius in Rhub
and Rtip
, respectively.
Assume the interfaces between elements are midway between adjacent element's midpoints.
AcousticAnalogies.source_elements_ccblade
— Methodsource_elements_ccblade(rotor::CCBlade.Rotor, sections::Vector{CCBlade.Section}, ops::Vector{CCBlade.OperatingPoint}, outputs::Vector{CCBlade.Outputs}, area_per_chord2::Vector{AbstractFloat}, period, num_src_times)
Construct and return an array of CompactSourceElement objects from CCBlade structs.
Arguments
rotor
: CCBlade rotor object.precone).sections
:Vector
of CCBlade section object.ops
:Vector
of CCBlade operating point.outputs
::Vector
of CCBlade output objects.area_per_chord2
: cross-sectional area divided by the chord squared of the element at each CCBlade.section. Should be a Vector{AbstractFloat}, same length assections
,ops
,outputs
.period
: length of the source time over which the returned source elements will evaluated.num_src_times
: number of source times.
AcousticAnalogies.to_paraview_collection
— Functionto_paraview_collection(name::AbstractString, ses::AbstractArray{<:CompactSourceElement}, axis::Integer=1)
Construct and write out a ParaView collection data file (.pvd
) object for an array of CompactSourceElement
with name name.pvd
(i.e., the name
argument should not contain a file extension).
axis
indicates the axis of ses
over which the source time for the source elements in ses
vary. One VTK PolyData (.vtp
) file will be written for each valid index along axis
.
Returns a list of filenames written out by WriteVTK.jl.
AcousticAnalogies.to_vtp
— Methodto_vtp(name::AbstractString, ses::AbstractArray{<:CompactSourceElement})
Construct and return a VTK polygonal (.vtp) data file object for an array of CompactSourceElement
with name name.vtp
(i.e., the name
argument should not contain a file extension).