General
DuctAPE.set_index_maps
— Functionset_index_maps(
npanels,
ncenterbody_inlet,
nwake_sheets,
dte_minus_cbte,
wnm,
wenids,
nwp,
nwsp,
nbn,
ndp,
riiw,
nrotor,
)
Set values for index map to be used throughout solve and post-process.
Arguments
npanels : paneling_constants.npanels
ncenterbody_inlet : paneling_constants.ncenterbody_inlet
nwake_sheets : paneling_constants.nwake_sheets
dte_minus_cbte : paneling_constants.dte_minus_cbte
wnm : wake_vortex_panels.nodemap
wenids : wake_vortex_panels.endnodeidxs
nwp : problem_dimensions.nwp
nwsp : problem_dimensions.nwsp
nbn : problem_dimensions.nbn
ndp : body_vortex_panels.npanel[1]
riiw : rotor_indices_in_wake
nrotor : problem_dimensions.nrotor
Returns
idmaps::NamedTuple
: A named tuple containing index mapping used in bookkeeping throughout solve and post-process
DuctAPE.precompute_parameters
— Functionprecompute_parameters(
ducted_rotor,
operating_point;
grid_solver_options=GridSolverOptions(),
integration_options=IntegrationOptions(),
autoshiftduct=true,
itcpshift=0.05,
axistol=1e-15,
tegaptol=1e1 * eps(),
finterp=(x,y,xp)->FLOWMath.akima(x,y,xp,2.0*eps(),eps()),
silence_warnings=true,
verbose=false,
)
Out of place main pre-processing function that computes all the required parameters for the solve.
Arguments
ducted_rotor::DuctedRotor
: A DuctedRotor objectoperating_point::OperatingPoint
: A OperatingPoint object
Keyword Arguments
grid_solver_options::GridSolverOptionsType=GridSolverOptions()
: A GridSolverOptionsType objectintegration_options::IntegrationMethod=IntegrationOptions()
: An IntegrationMethod objectautoshiftduct::Bool=true
: flag to shift duct geometry based on rotor tip radiusitcpshift::Float=0.05
: value used in positioning the internal pseudo control point in the solid bodies. Default is DFDC hard-coded value.axistol::Float=1e-15
: tolerance for how close to the axis of rotation to be considered on the axis.tegaptol::Float=1e1 * eps()
: tolerance for how large of a trailing edge gap is considered a gap.finterp::Function=FLOWMath.akima
: interpolation method for re-interpolating body coordinatessilence_warnings::Bool=true
: flag to silence warningsverbose::Bool=false
: flag to print verbose statements
Returns
ivr::NamedTuple
: A named tuple containing arrays of induced velocities on the rotorsivw::NamedTuple
: A named tuple containing arrays of induced velocities on the wakeivb::NamedTuple
: A named tuple containing arrays of induced velocities on the bodieslinsys::NamedTuple
: A named tuple containing cacheable data for the linear system, including:A_bb::Array{Float}
: AIC (LHS) matrix for the panel method systemb_bf::Array{Float}
: Initial system RHS vector based on freestrem magnitudeA_br::Array{Float}
: Unit normal velocity from rotors onto body panelsA_pr::Array{Float}
: Unit normal velocity from rotors onto body internal psuedo control pointsA_bw::Array{Float}
: Unit normal velocity from wake onto body panelsA_pw::Array{Float}
: Unit normal velocity from wake onto body internal psuedo control points
A_bb_LU::LinearAlgebra.LU
: LinearAlgebra LU factorization of the LHS matrixlu_decomp_flag::Vector{Bool}
: flag for whether factorization was successfulblade_elements::NamedTuple
: A named tuple containing cacheable blade element information (see docs forinterpolate_blade_elements
)airfoils::Vector{AFType}
: A matrix of airfoil types associated with each of the blade elementswakeK::Matrix{Float}
: A matrix of precomputed geometric constants used in the calculation of the wake vortex strengthsidmaps::NamedTuple
: A named tuple containing index mapping used in bookkeeping throughout solve and post-processpanels::NamedTuple
: A named tuple of panel objects including:body_vortex_panels::NamedTuple
: the named tuple containing the body vortex panel informationrotor_source_panels::NamedTuple
: the named tuple containing the rotor source panel informationwake_vortex_panels::NamedTuple
: the named tuple containing the wake vortex panel information
problem_dimensions::ProblemDimensions
: A ProblemDimensions object
precompute_parameters(
rp_duct_coordinates,
rp_centerbody_coordinates,
wake_grid,
rotor_indices_in_wake,
Rtips,
Rhubs,
rotor,
paneling_constants,
operating_point,
integration_options,
problem_dimensions=nothing;
itcpshift=0.05,
axistol=1e-15,
tegaptol=1e1 * eps(),
silence_warnings=true,
verbose=false,
)
An alternate version of precompute_parameters allowing for user defined geometry that does not go through a re-panling step (use with caution).
The first inputs are the outputs of the reinterpolate_geometry
and get_blade_ends_from_body_geometry
functions.
DuctAPE.precompute_parameters!
— Functionprecompute_parameters!(
ivr,
ivw,
blade_element_cache,
linsys,
wakeK,
ducted_rotor,
operating_point,
prepost_containers,
problem_dimensions;
grid_solver_options=GridSolverOptions(),
integration_options=IntegrationOptions(),
autoshiftduct=true,
itcpshift=0.05,
axistol=1e-15,
tegaptol=1e1 * eps(),
finterp=(x,y,xp)->FLOWMath.akima(x,y,xp,2.0*eps(),eps()),
silence_warnings=true,
verbose=false,
)
In-place version of precompute_parameters
.
precompute_parameters!(
ivr,
ivw,
blade_element_cache,
linsys,
wakeK,
wake_grid,
rp_duct_coordinates,
rp_centerbody_coordinates,
rotor_indices_in_wake,
rotor,
paneling_constants,
operating_point,
prepost_containers,
problem_dimensions=nothing;
integration_options=IntegrationOptions(),
itcpshift=0.05,
axistol=1e-15,
tegaptol=1e1 * eps(),
finterp=(x,y,xp)->FLOWMath.akima(x,y,xp,2.0*eps(),eps()),
silence_warnings=true,
verbose=false,
)
In-place version of the precompute_parameters
function by-passing the geometry reinterpolateion. (Use with caution)
Geometry
DuctAPE.reinterpolate_geometry
— Functionreinterpolate_geometry(
problem_dimensions,
duct_coordinates,
centerbody_coordinates,
rotor,
paneling_constants;
autoshiftduct=true,
grid_solver_options=GridSolverOptions(),
finterp=(x,y,xp)->FLOWMath.akima(x,y,xp,2.0*eps(),eps()),
verbose=false,
silence_warnings=true,
)
Re-interpolate the body geometry and return compatible body and way geometry.
Arguments
problem_dimensions::ProblemDimensions
: A ProblemDimensions objectduct_coordinates::Matrix{Float}
: [z,r] coordinates of duct geometrycenterbody_coordinates::Matrix{Float}
: [z,r] coordinates of centerbody geometryrotor::Rotor
: A Rotor objectpaneling_constants::PanelingConstants
: A PanelingConstants object
Keyword Arguments
autoshiftduct::Bool=true
: flag to shift duct geometry based on rotor tip radiusgrid_solver_options::SolverOptionsType=GridSolverOptions()
: options for the wake grid position solverfinterp::Function=FLOWMath.akima
: interpolation method for re-interpolating body coordinatesverbose::Bool=false
: flag to print verbose statementssilence_warnings::Bool=true
: flag to silence warnings
Returns
wake_grid::Array{Float}
: array containig the z and r elliptic grid points defning the wake geometry.rp_duct_coordinates::Matrix{Float}
: matrix containing the re-paneled duct coordinatesrp_centerbody_coordinates::Matrix{Float}
: matrix containing the re-paneled centerbody coordinatesrotor_indices_in_wake::Vector{Int}
: vector containing the indices of where in the wake the rotors reside (used later to define the rotor panel edges).
DuctAPE.reinterpolate_geometry!
— Functionreinterpolate_geometry!(
wake_grid,
rp_duct_coordinates,
rp_centerbody_coordinates,
rotor_indices_in_wake,
duct_coordinates,
centerbody_coordinates,
rotor,
blade_element_cache,
paneling_constants;
autoshiftduct=true,
grid_solver_options=GridSolverOptions(),
finterp=(x,y,xp)->FLOWMath.akima(x,y,xp,2.0*eps(),eps()),
verbose=false,
silence_warnings=true,
)
In-place version of reinterpolate_geometry
.
DuctAPE.generate_all_panels
— Functiongenerate_all_panels(
rp_duct_coordinates,
rp_centerbody_coordinates,
nwake_sheets,
rotor_indices_in_wake,
rotorzloc,
wake_grid;
itcpshift=0.05,
axistol=1e-15,
tegaptol=1e1 * eps(),
silence_warnings=true,
)
Function that calls all of the various panel generation functions are returns a named tuple containing all the panels
Arguments
rp_duct_coordinates::Matrix{Float}
: matrix containing the re-paneled duct coordinatesrp_centerbody_coordinates::Matrix{Float}
: matrix containing the re-paneled centerbody coordinatesnwake_sheets::Int
: number of wake sheetsrotor_indices_in_wake::Vector{Int}
: vector containing the indices of where in the wake the rotors reside (used later to define the rotor panel edges).rotorzloc:Vector{Float}
: axial locations of rotor lifting lines (contained in Rotor)wake_grid::Array{Float}
: array containig the z and r elliptic grid points defning the wake geometry.
Keyword Arguments
itcpshift::Float=0.05
: value used in positioning the internal pseudo control point in the solid bodies. Default is DFDC hard-coded value.axistol::Float=1e-15
: tolerance for how close to the axis of rotation to be considered on the axis.tegaptol::Float=1e1 * eps()
: tolerance for how large of a trailing edge gap is considered a gap.silence_warnings::Bool=true
: flag to silence warnings
Returns
panels::NamedTuple
: A named tuple of named tuples containing paneling information, including:body_vortex_panels::NamedTuple
rotor_source_panels::NamedTuple
wake_vortex_panels::NamedTuple
DuctAPE.generate_all_panels!
— Functiongenerate_all_panels!(
panels,
wake_grid,
rp_duct_coordinates,
rp_centerbody_coordinates,
rotor_indices_in_wake,
rotorzloc,
nwake_sheets;
itcpshift=0.05,
axistol=1e-15,
tegaptol=1e1 * eps(),
silence_warnings=true,
)
In-place version of generate_all_panels
.
Wake
DuctAPE.discretize_wake
— Functiondiscretize_wake(
duct_coordinates,
centerbody_coordinates,
rotorzloc, # rotor axial locations
wake_length,
npanels,
dte_minus_cbte;
)
Calculate wake sheet panel node z-coordinates.
Arguments
duct_coordinates::Matrix{Float}
: Array of input duct coordinatescenterbody_coordinates::Matrix{Float}
: Array of input centerbody_coordinates coordinatesrotorzloc ::Vector{Float}
: rotor axial locationswake_length::Float
: non-dimensional length of wake to extend beyond aft-most body trailing edge.npanels::Vector{Int}
: A vector of the number of panels between each discrete point. For example: [number of panels between the rotors; number of panels between the stator and the first trailing edge; number of panels between the trailing edges; number of panels between the last trailing edge and the end of the wake]dte_minus_cbte::Float
: indicator as to whether the duct trailing edge minus the centerbody trailing edge is positive, zero, or negative.
DuctAPE.generate_wake_grid
— Functiongenerate_wake_grid(
problem_dimensions,
rp_duct_coordinates,
rp_centerbody_coordinates,
Rhub1,
Rtip1,
tip_gap1,
zwake;
grid_solver_options=GridSolverOptions(),
verbose=false,
silence_warnings=true,
)
Initialize and solve for elliptic grid on which wake sheets are defined.
Arguments
problem_dimensions::
: A ProblemDimensions objectrp_duct_coordinates::
: repaneled duct coordinatesrp_centerbody_coordinates::
: repaneled centerbody coordinatesRhub1::
: Hub radius of first rotorRtip1::
: Tip radius of first rotortip_gap1::
: Tip gap of first rotor (MUST BE ZERO for now)zwake::
: axial positions of wake sheet panel nodes
Keyword Arguments
grid_solver_options::GridSolverOptionsType=GridSolverOptions()
: options for solving the elliptic grid.verbose::Bool=false
: flag to print verbose statementssilence_warnings::Bool=true
: flag to supress warnings
Returns
wake_grid::Array{Float,3}
: 3D Array of axial and radial wake_grid points after solution of elliptic system.
DuctAPE.generate_wake_grid!
— Functiongenerate_wake_grid!(
wake_grid,
rp_duct_coordinates,
rp_centerbody_coordinates,
Rhub1,
Rtip1,
tip_gap1,
zwake;
grid_solver_options=grid_solver_options,
verbose=false,
silence_warnings=true,
)
In-place version of generate_wake_grid
.
DuctAPE.initialize_wake_grid
— Functioninitialize_wake_grid(rp_duct_coordinates, rp_centerbody_coordinates, zwake, rwake)
Initialize the wake grid.
Arguments:
rp_duct_coordinates::Matrix{Float}
: The re-paneled duct coordinatesrp_centerbody_coordinates::Matrix{Float}
: The re-paneled centerbody coordinateszwake::Vector{Float}
: The axial positions of the wake sheet panel nodesrwake::Vector{Float}
: The radial positions of the blade elements for the foremost rotor
Returns:
wake_grid::Array{Float,3}
: 3D Array of axial and radial wake_grid points
DuctAPE.initialize_wake_grid!
— Functioninitialize_wake_grid!(
wake_grid, rp_duct_coordinates, rp_centerbody_coordinates, zwake, rwake
)
In-place version of initialize_wake_grid
.
DuctAPE.relax_grid!
— Functionrelax_grid!(
grid_solver_options::GridSolverOptionsType,
wake_grid;
verbose=false,
silence_warnings=true,
tabchar=" ",
ntab=1,
)
Relax/Solve initial wake grid according to elliptic system of equations.
Arguments
- `gridsolveroptions::GridSolverOptionsType' : options for elliptic grid solver
wake_grid::Array{Float,3}
: Initialized wake grid
Keyword Arguments
- `verbose=false::' : flag for printing verbose statements
- `silence_warnings=true::' : flag for supressing warnings
- `tabchar::String=" "::' : string to use for tabbing over verbose statements.
- `ntab::Int=1' : number of tabs for printing verbose statements
relax_grid!(xg, rg, nxi, neta; iteration_limit, atol)
Relax wakegrid using elliptic wakegrid solver.
Arguments:
wake_grid::Array{Float}
: initial guess for grid points.
Keyword Arguments:
iteration_limit::Int
: maximum number of iterations to run, default=100atol::Float
: convergence tolerance, default = 1e-9
DuctAPE.solve_elliptic_grid!
— Functionsolve_elliptic_grid!(
wake_grid;
algorithm=:trust_region,
autodiff=:forward,
atol=1e-14,
iteration_limit=10,
converged=[false],
residual_value=[0.0],
verbose=false,
)
Solve for elliptic grid using a non-SLOR approach that is compatible with ImplicitAD
Arguments:
wake_grid::Array{Float}
: initial guess for grid points.
Keyword Arguments:
algorithm::Symbol=:trust_region
: NLsolve algorithm typeautodiff::Symbol=:forward
: NLsolve derivatives optionatol::Float=1e-10
: convergence tolerance, default = 1e-9iteration_limit::Int=10
: maximum number of iterations to run, default=100converged::Vector{Bool}=[false]
: convergence flagresidual_value::Vector{Float}=[0.0]
: final residual valueverbose::Bool=false
: flag to print verbose statements
DuctAPE.solve_elliptic_grid
— Functionsolve_elliptic_grid(x,p)
Solve for elliptic grid using a non-SLOR approach that is compatible with ImplicitAD
Arguments:
x::Vector{Float}
: inputs (including initial guess for wake grid points)p::NamedTuple
: constant parameters
Returns:
converged_states::Vector{Float}
: converged states for grid solver
DuctAPE.elliptic_grid_residual!
— Functionelliptic_grid_residual!(r, y, x, p)
Residual function for elliptic grid solver.
Note: We assume grid has radial lines of constant axial position (as is done in the grid initialization) such that z'(ξ)=1.0 and z''(ξ)=z'(η)=z''(η)=0.0.
Arguments
r::Vector{Float}
: residual valuesy::Vector{Float}
: statesx::Vector{Float}
: inputs (to which derivatives are sensitive)p::Vector{Float}
: constants
DuctAPE.backward_stencil_1
— Functionbackward_stencil_1(f, x)
Three-point, backward difference scheme on a non-uniform grid.
Arguments:
f::AbstractMatrix{Float}
: f(x) for x{i-2}:x{i}x::AbstractVector{Float}
: [x{i-2}:x{i}]
Returns:
f'(x)::Float
: ∂f/∂x at x_i
DuctAPE.center_stencil_1
— Functioncenter_stencil_1(f, x)
Three-point, central difference scheme on a non-uniform grid, for first derivatives.
Arguments:
f::AbstractMatrix{Float}
: f(x) for x{i-1}:x{i+1}x::AbstractVector{Float}
: [x{i-1}:x{i+1}]
Returns:
f'(x)::Float
: ∂f/∂x at x_i
DuctAPE.center_stencil_2
— Functioncenter_stencil_2(f, x)
Three-point, central difference scheme on a non-uniform grid, for second derivatives.
Arguments:
f::AbstractMatrix{Float}
: f(x) for x{i-1}:x{i+1}x::AbstractVector{Float}
: [x{i-1}:x{i+1}]
Returns:
f''(x)::Float
: ∂^2f/∂x^2 at x_i
DuctAPE.center_stencil_2_mixed
— Functioncenter_stencil_2_mixed(f, x, y)
Three-point, central difference scheme on a non-uniform grid, for mixed derivatives. Uses a consecutive first derivative central difference (center_stencil_1
) in the x-dimension then another centeral difference of the f'(x) values in the y-dimension.
Arguments:
f::AbstractMatrix{Float}
: f(x,y) for x{i-1}:x{i+1} and y{i-1}:y{i+1}x::AbstractVector{Float}
: [x{i-1}:x{i+1}]y::AbstractVector{Float}
: [y{i-1}:y{i+1}]
Returns:
f'(x,y)::Float
: ∂^2f/∂x∂y at x_ij
DuctAPE.generate_wake_panels
— Functiongenerate_wake_panels(wake_grid)
Generate paneling for each wake sheet emanating from the rotor blade elements.
Arguments:
wake_grid::Array{Float,3}
: axial and radial locations of each wake_grid point (after relaxation/solution)
Returns:
wake_vortex_panels::NamedTuple
: A named tuple of panel values describing the wake vortex panels
DuctAPE.generate_wake_panels!
— Functiongenerate_wake_panels!(wake_panels, wake_grid)
In-place version of generate_wake_panels
.
DuctAPE.get_wake_k
— Functionget_wake_k(r, nwn)
Calculate geometric constant for use in later calculation of wake panel node strengths.
Arguments
r::Vector{Float}
: Vector of wake panel node radial positions
Returns
K::Vector{Float}
: Vector of geometric constants used in calculation of panel node strengths.
DuctAPE.get_wake_k!
— Functionget_wake_k!(K, r)
In-place version of get_wake_k
.
Bodies
DuctAPE.reinterpolate_bodies!
— Functionreinterpolate_bodies!(
rp_duct_coordinates,
rp_centerbody_coordinates,
duct_coordinates,
centerbody_coordinates,
zwake,
duct_le_coordinates,
ncenterbody_inlet,
nduct_inlet,
finterp=FLOWMath.akima,
)
Reinterpolate duct and centerbody coordinates in order to make them compatible with the calculated wake sheet panel axial positions.
Arguments
rp_duct_coordinates::Matrix{Float}
: the re-paneled duct coordinatesrp_centerbody_coordinates::Matrix{Float}
: the re-paneled centerbody coordinatesduct_coordinates::Matrix{Float}
: the input duct coordinatescenterbody_coordinates::Matrix{Float}
: the input centerbody coordinateszwake::Vector{Float}
: the wake sheet panel node axial positionsduct_le_coordinates::Matrix{Float}
: [z r] coordinates of duct leading edgencenterbody_inlet::Int
: the number of panels to use for the centerbody inletnduct_inlet,::Int
: the number of panels to use for the duct inlet
Keyword Arguments
finterp::Function=FLOWMath.akima
: interpolation method
DuctAPE.place_duct!
— Functionplace_duct!(rp_duct_coordinates, Rtip, rotorzloc, tip_gap)
Transform the duct radial coordinates such that the leading rotor radius touches the duct wall.
Note that this function is called AFTER the repanling function is called, such that the rotorzloc locations should line up directly with the duct and centerbody coordinates.
Arguments
rp_duct_coordinates::Matrix{Float}
: the re-paneled duct coordinatesRtip::Vector{Float}
: Tip radii for the rotor(s)rotorzloc::Vector{Float}
: axial position(s) of the rotor(s)tip_gap::Vector{Float}
: tip gap for the fore-most rotor (MUST BE ZERO for now)
Rotors
DuctAPE.interpolate_blade_elements
— Functioninterpolate_blade_elements(
rsp, Rtips, Rhubs, rotor_panel_centers, nbe; finterp=FLOWMath.linear
)
Interpolate blade elements based on Rotor inputs and number of desired blade elements (from number of wake sheet in PanelingConstants input)
Arguments
rsp::Rotor
: A Rotor object- `Rtips::Vector{Float}' : Vector of rotor tip radii
- `Rhubs::Vector{Float}' : Vector of rotor hub radii
- `rotorpanelcenters::Vector{Float}' : Vector of rotor panel centers
nbe::Int
: number of blade elements per rotor
Keyword Arguments
finterp::Function=FLOWMath.linear
: interpolation method (note, using Akima splines as is done for the body geometry can lead to negative chord in some cases)
Returns
blade_element_cache::NamedTuple
: A named tuple containing the cacheable blade element information excluding the airfoil data.airfoils::NamedTuple
: A named tuple containing vectors of inner and outer airfoil polar data for each blade element, used in interpolating the input data at blade element locations.
DuctAPE.interpolate_blade_elements!
— Functioninterpolate_blade_elements!(
blade_element_cache, rsp, rotor_panel_centers, nbe; finterp=FLOWMath.linear
)
In-place version of interpolate_blade_elements
.
Returns
airfoils::NamedTuple
: A named tuple containing vectors of inner and outer airfoil polar data for each blade element, used in interpolating the input data at blade element locations.
DuctAPE.get_blade_ends_from_body_geometry
— Functionget_blade_ends_from_body_geometry(
rp_duct_coordinates, rp_centerbody_coordinates, tip_gaps, rotorzloc
)
Obtain rotor hub and tip radii based on duct and centerbody geometry.
Arguments
var::type
:rp_duct_coordinates::Matrix{Float}
: re-paneled duct coordinatesrp_centerbody_coordinates::Matrix{Float}
: re-paneled centerbody coordinatestip_gaps::Vector{Float}
: gaps between blade tips and duct surface (MUST BE ZEROS for now)rotorzloc::Vector{Float}
: rotor lifting line axial positions.
Returns
Rtips::Vector{Float}
: rotor tip radiiRhubs::Vector{Float}
: rotor hub radii
DuctAPE.get_blade_ends_from_body_geometry!
— Functionget_blade_ends_from_body_geometry!(
Rtip,
Rhub,
rp_duct_coordinates,
rp_centerbody_coordinates,
tip_gaps,
rotorzloc;
silence_warnings=true,
)
In-place version of get_blade_ends_from_body_geometry
.
DuctAPE.get_local_solidity
— Functionget_local_solidity(B, chord, r)
Calculate local solidity from local chord, radial position, and number of blades.
Arguments
B::Float
: number of blades on rotor (usually an integer, but not necessarily).chord::Vector{Float}
: chord lengths at each radial station.r::Vector{Float}
: dimensional radial positions.
Returns
solidity::Vector{Float}
: local solidity at each radial station
DuctAPE.get_stagger
— Functionget_stagger(twists)
Convert twist angle to stagger angle
DuctAPE.generate_rotor_panels
— Functiongenerate_rotor_panels(rotorzloc, wake_grid, rotor_indices_in_wake, nwake_sheets)
Generate rotor panel objects.
Arguments
rotorzloc::Vector{Float}
: rotor lifting line axial positionwake_grid::Array{Float,3}
: wake elliptic grid axial and radial locationsrotor_indices_in_wake::Vector{Int}
: indices of where along wake the rotors are placednwake_sheets::Int
: number of wake sheets
Returns
rotor_source_panels::NamedTuple
: A named tuple containing the rotor source panel variables.
DuctAPE.generate_rotor_panels!
— Functiongenerate_rotor_panels!(
rotor_source_panels, rotorzloc, wake_grid, rotor_indices_in_wake, nwake_sheets
)
In-place version of generate_rotor_panels
.
Induced Velocities
DuctAPE.calculate_unit_induced_velocities
— Functioncalculate_unit_induced_velocities(problem_dimensions, panels, integration_options)
Calculate all the unit-induced velocties of all panels on all control points
Arguments
problem_dimensions::ProblemDimensions
: A ProblemDimensions objectpanels::NamedTuple
: A named tuple containing all the paneling informationintegration_options::IntegrationOptions
: Options used for integration of velocity kernals across panels
Returns
ivr::NamedTuple
: A named tuple containing arrays of induced velocities on the rotorsivw::NamedTuple
: A named tuple containing arrays of induced velocities on the wakeivb::NamedTuple
: A named tuple containing arrays of induced velocities on the bodies
DuctAPE.calculate_unit_induced_velocities!
— Functioncalculate_unit_induced_velocities!(ivr, ivw, ivb, panels, integration_options)
In-place version of calculate_unit_induced_velocities
.
DuctAPE.initialize_linear_system
— Functioninitialize_linear_system(
ivb,
body_vortex_panels,
rotor_source_panels,
wake_vortex_panels,
Vinf,
integration_options,
)
Set up the linear system used in the panel method solve.
Arguments
ivb::NamedTuple
: the named tuple containing all the unit induced velocities on the bodiesbody_vortex_panels::NamedTuple
: the named tuple containing the body vortex panel informationrotor_source_panels::NamedTuple
: the named tuple containing the rotor source panel informationwake_vortex_panels::NamedTuple
: the named tuple containing the wake vortex panel informationVinf::Vector{Float}
: the one-element vector containing the Freestream velocity magnitudeintegration_options::IntegrationOptions
: the integration options used in integrating the panel induced velocities
Returns
linsys::NamedTuple
: A named tuple containing cacheable data for the linear system, including:A_bb::Array{Float}
: AIC (LHS) matrix for the panel method systemb_bf::Array{Float}
: Initial system RHS vector based on freestrem magnitudeA_br::Array{Float}
: Unit normal velocity from rotors onto body panelsA_pr::Array{Float}
: Unit normal velocity from rotors onto body internal psuedo control pointsA_bw::Array{Float}
: Unit normal velocity from wake onto body panelsA_pw::Array{Float}
: Unit normal velocity from wake onto body internal psuedo control points
A_bb_LU::LinearAlgebra.LU
: LinearAlgebra LU factorization of the LHS matrixlu_decomp_flag::Vector{Bool}
: flag for whether factorization was successful
DuctAPE.initialize_linear_system!
— Functioninitialize_linear_system!(
linsys,
ivb,
body_vortex_panels,
rotor_source_panels,
wake_vortex_panels,
Vinf,
intermediate_containers,
integration_options,
)
In-place version of initialize_linear_system
.
Unit Induced Velocities
DuctAPE.calculate_xrm
— Functioncalculate_xrm(controlpoint, node)
Calculate xi, rho, and m for vortex and/or source ring induced velocity calculation.
Returns zeros if ring is on (or approximately on) the axis of rotation (zero radius).
Arguments
controlpoint::Vector{Float}
[z r] coordinates of point being influencednode::Vector{Float}
: [z r] coordinates of singularity ring
Returns
xi::Float
: normalized relative axial positionrho::Float
: normalized relative radial positionm::Float
: Elliptic integral inputrj::Float
: radial position of the ring
DuctAPE.calculate_xrm!
— Functioncalculate_xrm!(cache_vec, controlpoint, node)
In-place version of calculate_xrm
.
Cache_vec is a vector used to hold intermediate values as well as the outputs.
DuctAPE.get_elliptics
— Functionget_elliptics(m)
Calculate value of elliptic functions for the given geometry parameter.
Arguments
m::Float
: Elliptic Function parameter
Returns
K::Float
: K(m), value of elliptic function of the first kind at m.E::Float
: E(m), value of eeliptic function of the second kind at m.
DuctAPE.vortex_ring_vz
— Functionvortex_ring_vz(xi, rho, m, r_influence, influence_length)
Axial velocity induced by axisymmetric vortex ring.
Uses equivalent smoke ring induced velocity for self-induction, and returns zero if vortex ring is on axis of rotation (zero radius).
Arguments
xi::Float
: normalized z-coordinate, (z-zo)/rorho::Float
: normalized r-coordinate, r/rom::Float
: Elliptic Integral parameter, 4rho/sqrt(xi^2+(rho+1)^2)r_influence::Float
: radial location of vortex ring, roinfluence_length::Float
: length of panel used in calculating self-induction
Returns
vz::Float
: axially induced velocity of vortex ring
DuctAPE.vortex_ring_vz!
— Functionvortex_ring_vz!(xi, rho, m, r_influence, influence_length, cache_vec)
Same as vortexringvz, but uses the cache_vec to store intermediate calculations.
DuctAPE.smoke_ring_vz
— Functionsmoke_ring_vz(r_influence, influence_length)
Equivalent "smoke" ring self-induced velocity.
Arguments
r_influence::Float
: radial position of ring (i.e. the ring raidus)influence_length::Float
: length of influencing panel
Returs
vz::Float
: axially induced velocity of vortex ring
DuctAPE.vortex_ring_vr
— Functionvortex_ring_vr(xi, rho, m, r_influence)
Radial velocity induced by axisymmetric vortex ring.
Returns zero if vortex ring is on axis of rotation (zero radius), the point of influence is on the axis, or if self-inducing velocity.
Arguments
xi::Float
: normalized z-coordinate, (z-zo)/rorho::Float
: normalized r-coordinate, r/rom::Float
: Elliptic Integral parameter, 4rho/sqrt(xi^2+(rho+1)^2)r_influence::Float
: radial location of vortex ring, ro
Returns
vr::Float
: radially induced velocity of vortex ring
DuctAPE.vortex_ring_vr!
— Functionvortex_ring_vr!(xi, rho, m, r_influence, cache_vec)
Same as vortexringvr, but uses the cache_vec to store intermediate calculations.
DuctAPE.source_ring_vz
— Functionsource_ring_vz(xi, rho, m, r_influence)
Axial velocity induced by axisymmetric source ring.
Returns zero if source ring is on axis of rotation (zero radius), the point of influence is on the axis, or if self-inducing velocity.
Arguments:
xi::Float
: normalized z-coordinate, (z-zo)/rorho::Float
: normalized r-coordinate, r/rom::Float
: Elliptic Integral parameter, 4rho/sqrt(xi^2+(rho+1)^2)r_influence::Float
: radial location of vortex ring, ro
Returns:
vz::Float
: axially induced velocity of source ring
DuctAPE.source_ring_vz!
— Functionsource_ring_vz!(xi, rho, m, r_influence, cache_vec)
Same as sourceringvz, but uses cache_vec to store intermediate values.
DuctAPE.source_ring_vr
— Functionsource_ring_vr(xi, rho, m, r_influence)
Radial velocity induced by axisymmetric source ring.
Returns zero if source ring is on axis of rotation (zero radius), the point of influence is on the axis, or if self-inducing velocity.
Arguments:
xi::Float
: normalized z-coordinate, (z-zo)/rorho::Float
: normalized r-coordinate, r/rom::Float
: Elliptic Integral parameter, 4rho/sqrt(xi^2+(rho+1)^2)r_influence::Float
: radial location of vortex ring, ro
Returns:
vr::Float
: radially induced velocity of source ring
DuctAPE.source_ring_vr!
— Functionsource_ring_vr!(xi, rho, m, r_influence, cache_vec)
Same as sourceringvr, but uses cache_vec to store intermediate values.
Unit Induced Velocity Matrices
DuctAPE.induced_velocities_from_vortex_panels_on_points
— Functioninduced_velocities_from_vortex_panels_on_points(
controlpoints,
nodes,
nodemap,
influence_lengths,
integration_options;
integration_caches=nothing,
)
Calculate axial and radial components of induced velocity for a set of control points due to a set of axisymmetric vortex panels (bands).
Used for getting the unit induced velocities due to the body panels on the rotor/wake as well as the unit induced velocity due to the wake on the body/rotor.
Arguments
controlpoints::Matrix{Float}
[z r] coordinates of points being influencednodes::Matrix{Float}
: [z r] coordinates of vortex ringsnodemap::Matrix{Int}
: mapping from panel index to associated node indicesinfluence_lengths::Vector{Float}
: lengths over which vortex ring influence is applied on the surface.integration_options::IntegrationOptions
: integration options
Keyword Arguments
integration_caches::NamedTuple=nothing
: cache used in in-place quadrature functions.
Returns
VEL::Array{Float}
: N-controlpoint x N-node x [vz, vr] array of induced velocity components
DuctAPE.induced_velocities_from_vortex_panels_on_points!
— Functioninduced_velocities_from_vortex_panels_on_points!(
VEL,
controlpoint,
node,
nodemap,
influence_length,
integration_options;
integration_caches=nothing,
)
In-place version of induced_velocities_from_vortex_panels_on_points
.
DuctAPE.induced_velocities_from_source_panels_on_points
— Functioninduced_velocities_from_source_panels_on_points(
controlpoints,
nodes,
nodemap,
influence_lengths,
integration_options;
integration_caches=nothing,
)
Calculate axial and radial components of induced velocity for a set of control points due to a set of axisymmetric source panels (bands)
Used for getting the unit induced velocities due to the body panels on the rotor/wake as well as the unit induced velocity due to the wake on the body/rotor.
Arguments
controlpoints::Matrix{Float}
[z r] coordinates of points being influencednodes::Matrix{Float}
: [z r] coordinates of vortex ringsnodemap::Matrix{Int}
: mapping from panel index to associated node indicesinfluence_lengths::Vector{Float}
: lengths over which vortex ring influence is applied on the surface.integration_options::IntegrationOptions
: integration options
Returns:
VEL::Array{Float}
: N-controlpoint x N-node x [vz, vr] array of induced velocity components
DuctAPE.induced_velocities_from_source_panels_on_points!
— Functioninduced_velocities_from_source_panels_on_points!(
VEL,
controlpoint,
node,
nodemap,
influence_length,
integration_options;
integration_caches=nothing,
)
In-place version of induced_velocities_from_source_panels_on_points
.
DuctAPE.induced_velocities_from_trailing_edge_gap_panel!
— Functioninduced_velocities_from_trailing_edge_gap_panel!(
VEL,
controlpoint,
tenode,
teinfluence_length,
tendotn,
tencrossn,
teadjnodeidxs,
integration_options;
wake=false,
integration_caches=nothing,
)
Calculate axial and radial components of induced velocity for a set of control points due to any trailing edge gap panels.
Used for getting the unit induced velocities due to the body body trailing edge gap panels on the body/rotor/wake.
Note, this function is also used to calculate the influence of the wake ends rather than modeling a semi-infinite fortex sheet.
Arguments
VEL::Array{Float}
: N-controlpoint x N-node x [vz, vr] array of induced velocity components (modified in place)controlpoints::Matrix{Float}
[z r] coordinates of points being influencednodes::Matrix{Float}
: [z r] coordinates of vortex ringsnodemap::Matrix{Int}
: mapping from panel index to associated node indicesinfluence_lengths::Vector{Float}
: lengths over which vortex ring influence is applied on the surface.strengths::Matrix{Float}
: vortex constant circulation valuesintegration_options::IntegrationOptions
: integration options
Keyword Arguments
wake::Bool=false
: flag to indicate if this is being used for a wake sheet.integration_caches::NamedTuple=nothing
: cache used in in-place quadrature functions.
Panel Method Velocity Functions
DuctAPE.vortex_aic_boundary_on_boundary
— Functionvortex_aic_boundary_on_boundary(
controlpoint, normal, node, nodemap, influence_length, integration_options
)
Calculate panel method influence coefficients (V dot nhat) for a set of control points (on panels) due to a set of axisymmetric vortex rings (also on body surface)
Can be used for constructing the LHS influence Matrix for the panel method system.
Arguments
controlpoint::Matrix{Float}
[z r] coordinates of points being influencednormal::Matrix{Float}
: unit normal vectors of the panels on which the control points are situated.node::Matrix{Float}
: [z r] coordinates of panel nodes (edges)nodemap::Matrix{Int}
: [1 2] node indices for each panelinfluence_length::Vector{Float}
: lengths of influencing panelsintegration_options::IntegrationOptions
: integration options
Returns
AICn::Matrix{Float}
: N controlpoint x N+1 node array of V dot nhat values
DuctAPE.vortex_aic_boundary_on_boundary!
— Functionvortex_aic_boundary_on_boundary!(
AICn,
controlpoint,
normal,
node,
nodemap,
influence_length,
integration_options;
integration_caches=nothing,
)
In-place verion of vortex_aic_boundary_on_boundary
.
integration_caches is a named tuple containing caching for intermediate calculation values.
DuctAPE.vortex_aic_boundary_on_field
— Functionvortex_aic_boundary_on_field(
controlpoint,
normal,
node,
nodemap,
influence_length,
integration_options;
integration_caches=nothing,
)
Calculate panel method influence coefficients (V dot nhat) for a set of control points (NOT on panels) due to a set of axisymmetric vortex rings (on body surface)
Used for constructing portions of the panel method LHS matrix related to the pseudo control points in the bodies.
Arguments:
controlpoint::Matrix{Float}
[z r] coordinates of points being influencednormal::Matrix{Float}
: unit normal vectors of the panels on which the control points are situated.node::Matrix{Float}
: [z r] coordinates of panel nodes (edges)nodemap::Matrix{Int}
: [1 2] node indices for each panelinfluence_length::Vector{Float}
: lengths of influencing panelsintegration_options::IntegrationOptions
: integration options
Keyword Arguments
integration_caches::NamedTuple=nothing
: caches for intermediate values in integration.
Returns:
AICn::Matrix{Float}
: N controlpoint x N+1 node array of V dot nhat values
DuctAPE.vortex_aic_boundary_on_field!
— Functionvortex_aic_boundary_on_field!(
AICn,
controlpoint,
normal,
node,
nodemap,
influence_length,
integration_options;
integration_caches=nothing,
)
In-place version of vortex_aic_boundary_on_field
.
DuctAPE.add_kutta!
— Functionadd_kutta!(LHS, AICn, kids)
Add Kutta condition (γ1 + γN = 0) to LHS matrix.
LHS::Matrix{Float}
: a pre-allocated (zeros) full size left-hand side matrixAICn::Matrix{Float}
: influence coefficients for panels/nodeskids::Vector{Int}
: [1 2] indices of where to put 1's for kutta condition
DuctAPE.add_te_gap_aic!
— Functionadd_te_gap_aic!(
AICn,
controlpoint,
normal,
tenode,
teinfluence_length,
tendotn,
tencrossn,
teadjnodeidxs,
integration_options;
wake=false,
integration_caches=nothing,
)
Add trailing edge gap aerodynmic influence coefficient contributions to the AIC matrix.
Arguments
AICn::Matrix{Float}
: N controlpoint x N+1 node array of V dot nhat valuescontrolpoint::Matrix{Float}
[z r] coordinates of points being influencednormal::Matrix{Float}
: unit normal vectors of the panels on which the control points are situated.tenode::Matrix{Float}
: [z r] coordinates of trailing edge panel nodes (edges)teinfluence_length::Vector{Float}
: lengths of influencing trailing edge panelstendotn::Matrix{Float}
: nhat of trailing edge panel dotted with nhat of adjacent panelstencrossn::Matrix{Float}
: nhat of trailing edge panel crossed with nhat of adjacent panelsteadjnodeidxs::Matrix{Float}
: indices of nodes adjacent to trailing edge panelintegration_options::IntegrationOptions
: integration options
Keyword Arguments
wake::Bool=false
: flag as to whether this function is being applied to a wake sheet.integration_caches::NamedTuple=nothing
: caches for intermediate values in integration.
DuctAPE.source_aic
— Functionsource_aic(
controlpoint,
normal,
node,
nodemap,
influence_length,
integration_options;
integration_caches=nothing,
)
Calculate panel method influence coefficients (V dot nhat) for a set of control points (on panels) due to a set of axisymmetric source rings not on body surface.
Can be used for constructing the RHS boundary conditions due to rotor source panels.
Arguments
controlpoint::Matrix{Float}
[z r] coordinates of points being influencednormal::Matrix{Float}
: unit normal vectors of the panels on which the control points are situated.node::Matrix{Float}
: [z r] coordinates of panel nodes (edges)nodemap::Matrix{Int}
: [1 2] node indices for each panelinfluence_length::Vector{Float}
: lengths of influencing panelsintegration_options::IntegrationOptions
: integration options
Keyword Arguments
integration_caches::NamedTuple=nothing
: caches for intermediate values in integration.
Returns
AICn::Matrix{Float}
: N controlpoint x N+1 node array of V dot nhat values
DuctAPE.source_aic!
— Functionsource_aic!(
AICn,
controlpoint,
normal,
node,
nodemap,
influence_length,
integration_options;
integration_caches=nothing,
)
In-place version of source_aic
.
DuctAPE.freestream_influence_vector
— Functionfreestream_influence_vector(normals, Vinfmat)
Calculate RHS contributions due to freestream.
Note that the freestream is assumed to have zero radial component in the underlying theory, but here we allow an arbitrary 2D vector for velocity for taking the dot product easier.
Arguments
normals::Matrix{Float}
: unit normal vectors of the panels on which the control points are situated.Vinfmat::Matrix{Float}
: [z r] components of freestream velocity (r's should be zero)
Returns
RHS::Vector{Float}
: vector of normal components of freestream velocity on input panels
DuctAPE.freestream_influence_vector!
— Functionfreestream_influence_vector!(RHS, normals, Vinfmat)
In-place version of freestream_influence_vector
.
DuctAPE.assemble_lhs_matrix
— Functionassemble_lhs_matrix(
AICn, AICpcp, npanel, nnode, totpanel, totnode, prescribednodeidxs; dummyval=1.0
)
Assemble the LHS matrix of the panel method.
Arguments
AICn::Matrix{Float}
: N controlpoint x N+1 node array of V dot nhat valuesAICpcp::Matrix{Float}
: Nbodies controlpoint x N+1 node array of V dot nhat values (influence on psuedo control points)npanel::Vector{Int}
: number of panels comprising each bodynnode::Vector{Int}
: number of nodes comprising each bodytotpanel::Int
: total number of panelstotnode::Int
: total number of nodesprescribednodeidxs::Vector{Int}
: indices of nodes with prescribed strengths (those on the axis of rotation).
Keyword Arguments
dummyval::Float=1.0
: value for dummy input for prescribed and internal control points in the system. Do not change except for debugging purposes.
Returns
LHS::Matrix{Float}
: The full LHS matrix for the panel method.
DuctAPE.assemble_lhs_matrix!
— Functionassemble_lhs_matrix!(
LHS, AICn, AICpcp, npanel, nnode, totpanel, totnode, prescribednodeidxs; dummyval=1.0
)
In-place version of assemble_lhs_matrix
.
DuctAPE.factorize_LHS
— Functionfactorize_LHS(A::AbstractMatrix)
Returns the LU decomposition of A
.
DuctAPE.factorize_LHS!
— Functionfactorize_LHS!(Apivot::AbstractMatrix, A::AbstractMatrix)
Returns the LU decomposition of A
using Apivot
as storage memory to pivot leaving A
unchanged.
DuctAPE.assemble_rhs_matrix
— Functionassemble_rhs_matrix(
vdnb, vdnpcp, npanel, nnode, totpanel, totnode, prescribednodeidxs
)
Arguments
vdnb::Vector{Float}
: V dot nhat for body panelsvdnpcp::Vector{Float}
: V dot nhat for pseudo control pointsnpanel::Vector{Int}
: number of panels comprising each bodynnode::Vector{Int}
: number of nodes comprising each bodytotpanel::Int
: total number of body panelstotnode::Int
: total number of body nodesprescribednodeidxs::Vector{Int}
: indices of nodes with prescribed strengths (those on the axis of rotation)
Returns
RHS::Vector{Float}
: the RHS vector of the panel method.
DuctAPE.assemble_rhs_matrix!
— Functionassemble_rhs_matrix!(
RHS, vdnb, vdnpcp, npanel, nnode, totpanel, totnode, prescribednodeidxs
)
In-place version of assemble_rhs_matrix
.
DuctAPE.calculate_normal_velocity
— Functioncalculate_normal_velocity(velocity_vector, normal)
Calculate normal velocity_vector (V dot nhat).
Arguments
velocity_vector::Matrix{Float}
: velocity vector [z r] on each panelnormal::Matrix{Float}
: the panel unit normals
Returns
AIC::Matrix{Float}
: V dot n on each panel
DuctAPE.calculate_normal_velocity!
— Functioncalculate_normal_velocity!(AIC, velocity_vector, normal)
In-place version of calculate_normal_velocity
.
Quadrature
Integrands
DuctAPE.nominal_vortex_induced_velocity_sample
— Functionnominal_vortex_induced_velocity_sample(
t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0, 1.0)
)
Calculate the velocity induced by a linear vortex panel on a point.
Arguments
t::Float
: sample point in range (0,1) selected by quadrature method.node1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncache_vec::Vector{Float}
: cache for intermediate calculations
Keyword Arguments
nondimrange::Tuple=(0.0,1.0)
: Non-dimensional range describing the panel length. Do not change excpet for debugging purposes. Note, can also be a vector.
Returns
V::Matrix{Float}
: 2x2 matrix of axial and radial induced velocities from each of the nodes.
DuctAPE.nominal_vortex_induced_velocity_sample!
— Functionnominal_vortex_induced_velocity_sample!(
V, t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0, 1.0)
)
In-place version of nominal_vortex_induced_velocity_sample
.
DuctAPE.subtracted_singular_vortex_influence
— Functionsubtracted_singular_vortex_influence(node, controlpoint)
Calculate the singular portions of the self-induced vortex panel influence to subtract off the integral in the separation of singularity method.
Arguments
node::Vector{Float}
: node positioncontrolpoint::Vector{Float}
: controlpoint position
Returns
axial::Float
: axial direction influenceradial::Float
: radial direction influence
DuctAPE.subtracted_singular_vortex_influence!
— Functionsubtracted_singular_vortex_influence!(node, controlpoint, cache_vec)
Somewhat in-place version of subtracted_singular_vortex_influence
.
Arguments
node::Vector{Float}
: node positioncontrolpoint::Vector{Float}
: controlpoint positioncache_vec::Vector{Float}
: used to store intermediate values.
Returns
axial::Float
: axial direction influenceradial::Float
: radial direction influence
DuctAPE.analytically_integrated_vortex_influence
— Functionanalytically_integrated_vortex_influence(r, influence_length)
Analytical approximation of the singular portions of the self-induced vortex panel velocities to be added back in as part of the separation of singularity method.
Arguments
r::Float
: radial position of self-induced control pointinfluence_length::Float
: dimensional length of the panel
Returns
V::Vector{Float}
: axial and radial induced velocities
DuctAPE.analytically_integrated_vortex_influence!
— Functionanalytically_integrated_vortex_influence!(V, r, influence_length)
In-place version of analytically_integrated_vortex_influence
.
DuctAPE.self_vortex_induced_velocity_sample
— Functionself_vortex_induced_velocity_sample(
t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0, 1.0)
)
Calculate the velocity induced by a linear vortex panel on a point at the midpoint between the panel edges.
Arguments
t::Float
: sample point in range (0,1) selected by quadrature method.node1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncache_vec::Vector{Float}
: cache for intermediate calculations
Keyword Arguments
nondimrange::Tuple=(0.0,1.0)
: Non-dimensional range describing the panel length. Do not change excpet for debugging purposes. Note, can also be a vector.
Returns
V::Matrix{Float}
: 2x2 matrix of axial and radial induced velocities from each of the nodes.
DuctAPE.self_vortex_induced_velocity_sample!
— Functionself_vortex_induced_velocity_sample!(
V, t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0, 1.0)
)
In-place version of self_vortex_induced_velocity_sample
.
DuctAPE.nominal_source_induced_velocity_sample
— Functionnominal_source_induced_velocity_sample(
t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0, 1.0)
)
Calculate the velocity induced by a source panel on a point.
Arguments
t::Float
: sample point in range (0,1) selected by quadrature method.node1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncache_vec::Vector{Float}
: cache for intermediate calculations
Keyword Arguments
nondimrange::Tuple=(0.0,1.0)
: Non-dimensional range describing the panel length. Do not change excpet for debugging purposes. Note, can also be a vector.
Returns
V::Matrix{Float}
: 2x2 matrix of axial and radial induced velocities from each of the nodes.
DuctAPE.nominal_source_induced_velocity_sample!
— Functionnominal_source_induced_velocity_sample!(
V, t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0,1.0)
)
In-place version of nominal_source_induced_velocity_sample
.
DuctAPE.subtracted_singular_source_influence
— Functionsubtracted_singular_source_influence(node, controlpoint)
Calculate the singular portions of the self-induced source panel influence to subtract off the integral in the separation of singularity method.
Arguments
node::Vector{Float}
: node positioncontrolpoint::Vector{Float}
: controlpoint position
Returns
axial::Float
: axial direction influenceradial::Float
: radial direction influence
DuctAPE.subtracted_singular_source_influence!
— Functionsubtracted_singular_source_influence!(node, controlpoint, cache_vec)
In-place version of subtracted_singular_source_influence
.
DuctAPE.analytically_integrated_source_influence
— Functionanalytically_integrated_source_influence(r, influence_length)
In-place version of analytically_integrated_source_influence
.
DuctAPE.analytically_integrated_source_influence!
— Functionanalytically_integrated_source_influence(r, influence_length)
Analytical approximation of the singular portions of the self-induced source panel velocities to be added back in as part of the separation of singularity method.
Arguments
r::Float
: radial position of self-induced control pointinfluence_length::Float
: dimensional length of the panel
Returns
V::Vector{Float}
: axial and radial induced velocities
DuctAPE.self_source_induced_velocity_sample
— Functionself_source_induced_velocity_sample(
t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0, 1.0)
)
Calculate the velocity induced by a linear source panel on a point at the midpoint between the panel edges.
Arguments
t::Float
: sample point in range (0,1) selected by quadrature method.node1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncache_vec::Vector{Float}
: cache for intermediate calculations
Keyword Arguments
nondimrange::Tuple=(0.0,1.0)
: Non-dimensional range describing the panel length. Do not change excpet for debugging purposes. Note, can also be a vector.
Returns
V::Matrix{Float}
: 2x2 matrix of axial and radial induced velocities from each of the nodes.
DuctAPE.self_source_induced_velocity_sample!
— Functionself_source_induced_velocity_sample!(
V, t, node1, node2, influence_length, controlpoint, cache_vec; nondimrange=(0.0, 1.0)
)
In-place version of self_source_induced_velocity_sample
.
Integrals
DuctAPE.nominal_vortex_panel_integration
— Functionnominal_vortex_panel_integration(
integration_options,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
Integration of vortex panel induced velocity on a point far away.
Arguments
integration_options::IntegrationMethod
: options for itegration methodsnode1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncontainers::NamedTuple
: cache for intermediate calculations
Keyword Arguments
debug::Bool=false
: if true, some methods will return the estimation error.
Returns
V::Matrix{Float}
: velocity components due to the jth and j+1th nodes in the format:[vz_j vr_j; vz_{j+1} vr_{j+1}]
DuctAPE.nominal_vortex_panel_integration!
— Functionnominal_vortex_panel_integration!(
integration_options,
V,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
In-place version of nominal_vortex_panel_integration
.
DuctAPE.self_vortex_panel_integration
— Functionself_vortex_panel_integration(
integration_options,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
Integration of linear vortex panel self-induced velocity.
Arguments
integration_options::IntegrationMethod
: options for itegration methodsnode1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncontainers::NamedTuple
: cache for intermediate calculations
Keyword Arguments
debug::Bool=false
: if true, some methods will return the estimation error.
Returns
V::Matrix{Float}
: velocity components due to the jth and j+1th nodes in the format:[vz_j vr_j; vz_{j+1} vr_{j+1}]
DuctAPE.self_vortex_panel_integration!
— Functionself_vortex_panel_integration!(
integration_options,
V,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
In-place version of self_vortex_panel_integration
.
DuctAPE.nominal_source_panel_integration
— Functionnominal_source_panel_integration(
integration_options,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
Integration of source panel induced velocity on a point far away.
Arguments
integration_options::IntegrationMethod
: options for itegration methodsnode1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncontainers::NamedTuple
: cache for intermediate calculations
Keyword Arguments
debug::Bool=false
: if true, some methods will return the estimation error.
Returns
V::Matrix{Float}
: velocity components due to the jth and j+1th nodes in the format:[vz_j vr_j; vz_{j+1} vr_{j+1}]
DuctAPE.nominal_source_panel_integration!
— Functionnominal_source_panel_integration!(
integration_options,
V,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
In-place version of nominal_source_panel_integration
.
DuctAPE.self_source_panel_integration
— Functionself_source_panel_integration(
integration_options,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
Integration of linear source panel self-induced velocity.
Arguments
integration_options::IntegrationMethod
: options for itegration methodsnode1::Vector{Float}
: first panel node (edge) position.node2::Vector{Float}
: second panel node (edge) position.influence_length::Float
: dimensional length of panel.controlpoint::Vector{Float}
: controlpoint positioncontainers::NamedTuple
: cache for intermediate calculations
Keyword Arguments
debug::Bool=false
: if true, some methods will return the estimation error.
Returns
V::Matrix{Float}
: velocity components due to the jth and j+1th nodes in the format:[vz_j vr_j; vz_{j+1} vr_{j+1}]
DuctAPE.self_source_panel_integration!
— Functionself_source_panel_integration!(
integration_options,
V,
node1,
node2,
influence_length,
controlpoint,
containers;
debug=false,
)
In-place version of self_source_panel_integration
.
DuctAPE.extrapolate!
— Functionextrapolate!(V, err, fh; power=2, atol=1e-6)
Performs Richardson extrapolation on an array fh
for use in Romberg integration.
Arguments
V::Matrix{Float}
: velocity components due to the jth and j+1th nodes in the format:[vz_j vr_j; vz_{j+1} vr_{j+1}]
err::Vector{Float}
: estimated errors in velocity approximationfh::Tuple
:(f(h), h)
tuples (in order of decreasing|h|
)
State Initialization
DuctAPE.initialize_velocities
— Functioninitialize_velocities(
solver_options::SolverOptionsType,
operating_point,
blade_elements,
linsys,
ivr,
ivw,
body_totnodes,
wake_panel_sheet_be_map,
)
Initialize velocity state variables.
Arguments
solver_options::SolverOptionsType
: solver options type for dispatchoperating_point::OperatingPoint
: an OperatingPoint objectblade_elements::NamedTuple
: A named tuple containing the blade element geometry and airfoil information.linsys::NamedTuple
: A named tuple containing the panel method linear system information.ivr::NamedTuple
: A named tuple containing the unit induced velocities on the rotorsivw::NamedTuple
: A named tuple containing the unit induced velocities on the wakebody_totnodes::Int
: the total number of panel nodes comprising the duct and centerbody geometrywake_panel_sheet_be_map::Matrix{Int}
: An index map from the wake panels to the nearest ahead rotor blade element along the wake sheets
Returns
vz_rotor::Vector{Float}
: a vector of the velocity state variables associated with the rotor axially induced velocityvtheta_rotor::Vector{Float}
: a vector of the velocity state variables associated with the rotor tangentially induced velocityCm_wake::Vector{Float}
: a vector of the velocity state variables associated with the wake control point meridional velocity
DuctAPE.initialize_velocities!
— Functionfunction initializevelocities!( solveroptions::SolverOptionsType, vzrotor, vthetarotor, Cmwake, operatingpoint, bladeelements, linsys, ivr, ivw, bodytotnodes, wakepanelsheetbemap, )
In-place version of initialize_velocities
.
DuctAPE.initialize_strengths!
— Functioninitialize_strengths!(
solver_options::SolverOptionsType,
Gamr,
sigr,
gamw,
operating_point,
blade_elements,
linsys,
ivr,
ivw,
wakeK,
body_totnodes,
wake_nodemap,
wake_endnodeidxs,
wake_panel_sheet_be_map,
wake_node_sheet_be_map,
wake_node_ids_along_casing_wake_interface,
wake_node_ids_along_centerbody_wake_interface,
)
Initialize strength state variables.
Arguments
solver_options::SolverOptionsType
: solver options type for dispatchGamr::Vector{Float}
: Rotor circulation state variables (modified in place)sigr::Vector{Float}
: Rotor panel strength state variables (modified in place)gamw::Vector{Float}
: Wake panel strength state variables (modified in place)operating_point::OperatingPoint
: an OperatingPoint objectblade_elements::NamedTuple
: A named tuple containing the blade element geometry and airfoil information.linsys::NamedTuple
: A named tuple containing the panel method linear system information.ivr::NamedTuple
: A named tuple containing the unit induced velocities on the rotorsivw::NamedTuple
: A named tuple containing the unit induced velocities on the wakewakeK::Vector{Float}
: geometric constants of wake nodes used in calculating wake strengthsbody_totnodes::Int
: the total number of panel nodes comprising the duct and centerbody geometrywake_nodemap::Matrix{Int}
: an index map of wake panel to the associated node indiceswake_endnodeidxs::Matrix{Int}
: the node indices of the start and end of the wake sheets.wake_panel_sheet_be_map::Matrix{Int}
: An index map from the wake panels to the nearest ahead rotor blade element along the wake sheetswake_node_sheet_be_map::Matrix{Int}
: An index map from the wake nodes to the nearest ahead rotor blade element along the wake sheetswake_node_ids_along_casing_wake_interface::type
: An index map indicating which wake nodes interface with the duct wallwake_node_ids_along_centerbody_wake_interface::type
: An index map indicating which wake nodes interface with the centerbody wall
function initialize_strengths!(
solver_options::CSORSolverOptions,
Gamr,
sigr,
gamw,
solve_containers,
operating_point,
blade_elements,
wakeK,
wake_nodemap,
wake_endnodeidxs,
wake_panel_sheet_be_map,
wake_node_sheet_be_map,
wake_node_ids_along_casing_wake_interface,
wake_node_ids_along_centerbody_wake_interface;
niter=10,
rlx=0.5,
)
Refactored from DFDC SUBROUTINE ROTINITBLD
From the subroutine notes: Sets reasonable initial circulation using current rotor blade geometry (chord, beta). Initial circulations are set w/o induced effects An iteration is done using the self-induced velocity from momentum theory to converge an approximate induced axial velocity
Arguments
solver_options::SolverOptionsType
: solver options type for dispatchGamr::Vector{Float}
: Rotor circulation state variables (modified in place)sigr::Vector{Float}
: Rotor panel strength state variables (modified in place)gamw::Vector{Float}
: Wake panel strength state variables (modified in place)operating_point::OperatingPoint
: an OperatingPoint objectblade_elements::NamedTuple
: A named tuple containing the blade element geometry and airfoil information.wakeK::Vector{Float}
: geometric constants of wake nodes used in calculating wake strengthswake_nodemap::Matrix{Int}
: an index map of wake panel to the associated node indiceswake_endnodeidxs::Matrix{Int}
: the node indices of the start and end of the wake sheets.wake_panel_sheet_be_map::Matrix{Int}
: An index map from the wake panels to the nearest ahead rotor blade element along the wake sheetswake_node_sheet_be_map::Matrix{Int}
: An index map from the wake nodes to the nearest ahead rotor blade element along the wake sheetswake_node_ids_along_casing_wake_interface::type
: An index map indicating which wake nodes interface with the duct wallwake_node_ids_along_centerbody_wake_interface::type
: An index map indicating which wake nodes interface with the centerbody wall
Keyword Arguments
rlx::Float=0.5
: factor for under-relaxation to reduce transients in CL
Returns