Public API

Input Types

DuctAPE.PanelingConstantsType
PanelingConstants(
    nduct_inlet,
    ncenterbody_inlet,
    npanels,
    dte_minus_cbte,
    nwake_sheets,
    wake_length=1.0,
)

Constants used in re-paneling geometry.

Note that unlike other input structures, this one, in general, does not define fields as vectors. This is because these values should not change throughout an optimization, even if the geometry may change. Otherwise, discontinuities could be experienced.

Arguments

  • nduct_inlet::Int : The number of panels to use for the casing inlet.
  • ncenterbody_inlet::Int : The number of panels to use for the centerbody inlet.
  • npanels::AbstractVector{Int} : A vector containing the number of panels between discrete locations inside the wake. Specifically, the number of panels between the rotors, between the last rotor and the first body trailing edge, between the body trailing edges (if different), and between the last body trailing edge and the end of the wake. The length of this vector should be N+1 (where N is the number of rotors) if the duct and centerbody trailing edges are aligned, and N+2 if not.
  • dte_minus_cbte::Float : An indicator concerning the hub and duct trailing edge relative locations. Should be set to -1 if the duct trailing edge axial position minus the centerbody trailing edge axial position is negative, +1 if positive (though any positive or negative number will suffice), and zero if the trailing edges are aligned.
  • nwake_sheets::Int : The number of wake sheets to use. Note this will also be setting the number of blade elements to use.
  • wake_length::Float=1.0 : Non-dimensional (based on the length from the foremost body leading edge and the aftmost body trailing edge) length of the wake extending behind the aftmost body trailing edge.
source
DuctAPE.RotorType
Rotor(
    B, rotorzloc, r, Rhub, Rtip, chords, twists, tip_gap, airfoils, fliplift
)

Composite type containing the rotor(s) geometric properties.

Note that the actual struct requires the inputs to be arrays, but there is a constructor function that will take in scalars and automatically build the array-based struct.

Arguments

  • B::AbstractVector{Float} : The number of blades for each rotor. May not be an integer, but usually is.
  • rotorzloc::AbstractVector{Float} : Dimensional, axial position of each rotor.
  • r::AbstractArray{Float} : Non-dimensional radial locations of each blade element.
  • Rhub::AbstractVector{Float} : Dimensional hub radius of rotor. (may be changed if it does not match the radial position of the centerbody geometry at the selected rotorzloc.
  • Rtip::AbstractVector{Float} : Dimensional tip radius of rotor. Is used to determine the radial position of the duct if the autoshiftduct option is selected.
  • chords::AbstractArray{Float} : Dimensional chord lengths of the blade elements.
  • twists::AbstractArray{Float} : Blade element angles, in radians.
  • tip_gap::AbstractVector{Float} : Currently unused, do not set to anything other than zeros.
  • airfoils::AbstractArray{AFType} : Airfoil types describing the airfoil polars for each blade element. Currently only fully tested with C4Blade.DFDCairfoil types.
  • fliplift::AbstractVector{Bool} : Flag to indicate if the airfoil lift values should be flipped or not.
source
DuctAPE.DuctedRotorType
DuctedRotor(duct_coordinates, centerbody_coordinates, rotor, paneling_constants)

Arguments

  • duct_coordinates::AbstractMatrix : The [z, r] coordinates of the duct geometry beginning at the inner (casing) side trailing edge and proceeding clockwise. Note that the duct geometry absolute radial position does not need to be included here if the autoshiftduct option is selected.
  • centerbody_coordinates::AbstractMatrix : The [z, r] coordinates of the centerbody beginning at the leading edge and ending at the trailing edge. Note that the leading edge is assumed to be placed at a radial distance of 0.0 from the axis of rotation.
  • rotor::Rotor : Rotor (and possibly stator) geometric paramters.
  • paneling_constants::PanelingConstants : Constants used in re-paneling the geometry.
source
DuctAPE.OperatingPointType
OperatingPoint(Vinf, Minf, rhoinf, muinf, asound, Ptot, Ttot, Omega)
OperatingPoint(
    Vinf, Omega, rhoinf=nothing, muinf=nothing, asound=nothing; altitude=0.0
)
OperatingPoint(
    ::Imperial, Vinf, Omega, rhoinf=nothing, muinf=nothing, asound=nothing; altitude=0.0
)

DuctedRotor operating point information.

Functions that take in altitude will populate undefined thermodynamic properties of the freestream using a standard_atmosphere model, ideal gas law, and Sutherland's law; defaulting to SI units. If the ::Imperial dispatch type is input, then the thermodynamic properties will be converted to Imperial units.

Fields/Arguments

  • Vinf::AbstractVector{Float} : Freestream velocity magnitude (which is only in the axial direction).
  • Minf::AbstractVector{Float} : Freestream Mach number
  • rhoinf::AbstractVector{Float} : Freestream density
  • muinf::AbstractVector{Float} : Freestream viscosity
  • asound::AbstractVector{Float} : Freestream speed of sound
  • Ptot::AbstractVector{Float} : Freestream total pressure
  • Ttot::AbstractVector{Float} : Freestream total temperature
  • Omega::AbstractVector{Float} : Rotor rototation rate(s)
source
DuctAPE.ReferenceParametersType
ReferenceParameters(Vref, Rref)

Reference parameters for post-process non-dimensionalization.

Note that the actual struct requires the inputs to be arrays, but there is a constructor function that will take in scalars and automatically build the array-based struct.

Arguments

  • Vref::AbstractVector{Float} : Reference velocity.
  • Rref::AbstractVector{Float} : Reference rotor tip radius.
source

Preallocations

DuctAPE.allocate_prepost_container_cacheFunction
allocate_prepost_container_cache(paneling_constants::PanelingConstants)
allocate_prepost_container_cache(problem_dimensions::ProblemDimensions)

Allocate the pre- and post-processing cache (used for intermediate calculations) based on paneling constants or problem dimensions.

Arguments

  • paneling_constants::PanelingConstants : a PanelingConstants object

OR

  • problem_dimensions::ProblemDimensions : a ProblemDimensions object

Keyword Arguments

  • fd_chunk_size::Int=12 : chunk size to use for PreallocationTools caches. Note that the automated chunk size for DuctAPE will always be the ForwardDiff threshold of 12 due to the size of the system, so it will be best to leave this at the default unless further development allows for chunk size selection for individual solvers.
  • levels::Int=1 : levels for nested duals. Note that since ImplicitAD is being used for all solves, there should be no need for more than 1 level.

Returns

  • prepost_container_caching::NamedTuple : a Named Tuple containing:
    • prepost_container_cache::PreallocationTools.DiffCache : the cache
    • prepost_container_cache_dims::NamedTuple : a named tuple containing the dimensions used for reshaping the cache when needed.
source
DuctAPE.allocate_solve_parameter_cacheFunction
allocate_solve_parameter_cache(
    solve_type::SolverOptionsType,
    paneling_constants::PanelingConstants;
    fd_chunk_size=12,
    levels=1,
)
allocate_solve_parameter_cache(
    solve_type::SolverOptionsType,
    problem_dimensions::ProblemDimensions;
    fd_chunk_size=12,
    levels=1
)

Allocate the solve parameter cache for parameters passed into the solver(s).

Arguments

  • solve_type::SolverOptionsType : Solver options type used for dispatch
  • paneling_constants::PanelingConstants : a PanlingConstants object used for sizing

OR

  • problem_dimensions::ProblemDimensions : a ProblemDimensions object used for sizing

Keyword Arguments

  • fd_chunk_size::Int=12 : chunk size to use for PreallocationTools caches. Note that the automated chunk size for DuctAPE will always be the ForwardDiff threshold of 12 due to the size of the system, so it will be best to leave this at the default unless further development allows for chunk size selection for individual solvers.
  • levels::Int=1 : levels for nested duals. Note that since ImplicitAD is being used for all solves, there should be no need for more than 1 level.

Returns

  • solve_parameter_caching::NamedTuple : a Named Tuple containing:
    • solve_parameter_cache::PreallocationTools.DiffCache : the cache
    • solve_parameter_cache_dims::NamedTuple : a named tuple containing the dimensions used for reshaping the cache when needed.
source
DuctAPE.allocate_solve_container_cacheFunction
allocate_solve_container_cache(
    solve_type::SolverOptionsType,
    paneling_constants::PanelingConstants;
    fd_chunk_size=12,
    levels=1,
)
allocate_solve_container_cache(
    solve_type::SolverOptionsType,
    problem_dimensions::ProblemDimensions;
    fd_chunk_size=12,
    levels=1,
)

Allocate the solve cache (used for intermediate calculations) based on paneling constants or problem dimensions.

Arguments

  • paneling_constants::PanelingConstants : a PanelingConstants object

OR

  • problem_dimensions::ProblemDimensions : a ProblemDimensions object

Keyword Arguments

  • fd_chunk_size::Int=12 : chunk size to use for PreallocationTools caches. Note that the automated chunk size for DuctAPE will always be the ForwardDiff threshold of 12 due to the size of the system, so it will be best to leave this at the default unless further development allows for chunk size selection for individual solvers.
  • levels::Int=1 : levels for nested duals. Note that since ImplicitAD is being used for all solves, there should be no need for more than 1 level.

Returns

  • solve_container_caching::NamedTuple : a Named Tuple containing:
    • solve_container_cache::PreallocationTools.DiffCache : the cache
    • solve_container_cache_dims::NamedTuple : a named tuple containing the dimensions used for reshaping the cache when needed.
source

Options

General Options

DuctAPE.OptionsType
struct Options

Type containing (nearly) all the available user options.

Fields

General Options

  • verbose::Bool = false : flag to print verbose statements
  • silence_warnings::Bool = true : flag to silence warnings
  • multipoint_index::Int = [1] : holds current index of multi-point solver (no need for user to change this usually)

Pre-processing Options

Geometry interpolation and generation options :

  • finterp::Interplation Method = FLOWMath.akima : interpolation method used for re-paneling bodies
  • autoshiftduct::Bool = true : flag as to whether duct geometry should be shifted based on rotor tip location
  • lu_decomp_flag::Bool = false : flag indicating if panel method LHS matrix factorization was successful

paneling options

  • itcpshift::Float = 0.05 : factor for internal trailing edge psuedo-panel placement (default is DFDC hard-coded value)
  • axistol::Float = 1e-15 : tolerance for how close the the axis of rotation should be considered on the axis
  • tegaptol::Float = 1e1 * eps() : tolerance for how large of a trailing edge gap should be considered a gap

Integration Options

  • integration_options::IntegrationOptions type = IntegrationOptions() : integration options

Post-processing Options

  • boundary_layer_options::BoundaryLayerOptions : BoundaryLayerOptions object
  • write_outputs::AbstractArray{Bool} = [false] : Bool for whether to write the outputs of the analysis to an external file (slow)
  • outfile::AbstractArray{String} = ["outputs.jl"] : External output file name (including path information) for files to write
  • checkoutfileexists::Bool = false : Flag for whether to check if file exists before overwriting
  • output_tuple_name::AbstractArray{String} = ["outs"] : variable name for named tuple written to out file

Solving Options

  • grid_solver_options::GridSolverOptionsType = GridSolverOptions() : elliptic grid solver options
  • solver_options::SolverOptionsType = ChainSolverOptions() : solver options

Failure Options

  • hard_fail::Bool = true : flag as to whether DuctAPE should return nothing immediately after a failed initialization of the elliptic grid or a failed decomposition of the body influence matrix. If set to false, DuctAPE will attempt to return objects of the correct size, but with initialized values only.
source
DuctAPE.set_optionsFunction
set_options(; kwargs...)
set_options(multipoint; kwargs...)

Set the options for DuctAPE to use.

Note that the vast majority of the available options are defined through keyword arguments. See the documentation for the various option types for more information.

Arguments

  • multipoint::AbstractArray{OperatingPoint} : a vector of operating points to use if running a multi-point analysis.
source

Integration Options

DuctAPE.IntegrationOptionsType
struct IntegrationOptions

A struct used to hold the integration options for both the nominal and singular cases.

Fields

  • nominal::IntegrationMethod=GaussLegendre(8) : the integration options to use for the nominal case.
  • singular::IntegrationMethod=GaussLegendre(8) : the integration options to use for the self-induced case.
source
DuctAPE.GaussLegendreType
struct GaussLegendre <: IntegrationMethod

Options for Gauss-Legendre integration method

Fields

  • sample_points::Int : Sample Points
  • weights::Int : Gauss weights
source
DuctAPE.GaussKronrodType
struct GaussKronrod <: IntegrationMethod

Options for Gauss-Kronrod integration method

Fields

  • order::Int = 7 : order of Legendre polynomial to use on each interval
  • maxevales::Int = 10^7 : maximum number of evaluations in the adaptive method
  • atol::Float = 0.0 : absolute error tolerance. (note, if zero, QuadGK uses sqrt(eps()) relative tolerance).
source
DuctAPE.RombergType
struct Romberg <: IntegrationMethod

Options for Romberg integration method

Fields

  • max_subdivisions::Int = 10 : maximum number of subdivisions. Note, total number of internvals is 2^N, where N is number of subdivisions.
  • atol::Float = 1e-6 : absolute error tolerance.
source

Solver Options

Elliptic Grid Solve

DuctAPE.SLORGridSolverOptionsType
struct SLORGridSolverOptions <: GridSolverOptionsType

Options for SLOR (successive line over relaxation) elliptic grid solver.

Fields

  • iteration_limit::Int = 100 : maximum number of iterations
  • atol::Float = 1e-9 : absolute convergence tolerance
  • converged::AbstractArray{Bool} = [false]
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.GridSolverOptionsType
struct GridSolverOptions <: GridSolverOptionsType

Options for Newton elliptic grid solver.

Fields

  • iteration_limit::Int = 20 : maximum number of iterations
  • atol::Float = 3e-10 : absolute convergence tolerance
  • algorithm::Symbol = :newton : algorithm to use in NLsolve.jl
  • autodiff::Symbol = :forward : differentiation method to use in NLsolve.jl
  • precondition = false : flag to precondition with SLOR
  • precondition_max_iterations = 3 : number of precondition iterations
  • converged::AbstractArray{Bool} = [false]
  • iterations::AbstractArray{Int} = [0] : iteration counter
  • residual_value::AbstractArray{Int} = [0] : residual value
source

Aerodynamics Solve

DuctAPE.ChainSolverOptionsType
struct ChainSolverOptions <:ExternalPolyAlgorithmOptions

Options for Chain Solvers (try one solver, if it doesn't converge, try another)

Fields

  • `solvers::AbstractArray{SolverOptionsType} = [ NLsolveOptions(; algorithm=:anderson, atol=1e-12), MinpackOptions(; atol=1e-12), NonlinearSolveOptions(; algorithm=SimpleNonlinearSolve.SimpleNewtonRaphson, atol=1e-12, additional_kwargs=(; autodiff=SimpleNonlinearSolve.AutoForwardDiff()), ), ] : Vector of solver options to use.
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.CompositeSolverOptionsType
struct CompositeSolverOptions <: ExternalPolyAlgorithmOptions

Options for Composite Solvers (start with a partial solve of one solve, then finish with another starting where the first left off).

Fields

  • `solvers::AbstractArray{SolverOptionsType} = [ NLsolveOptions(; algorithm=:newton, iteration_limit=3), NLsolveOptions(; algorithm=:anderson, atol=1e-12), ]' : Vector of solver options to use.
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.NLsolveOptionsType
struct NLsolveOptions <: ExternalSolverOptions

Options for the NLsolve pacakge solvers

Fields

  • algorithm::Symbol = :anderson : algorithm to use
  • additional_kwargs = (;) : any additional keyword arguments for the solver
  • atol::Float = 1e-12 : absolute convergence tolerance
  • iteration_limit::Int = 25 : maximum number of iterations
  • linesearch_method::LineSearches method = LineSearches.MoreThuente : line search method to use
  • linesearch_kwargs = (;) : any additional lineseach keyword arguments
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.NonlinearSolveOptionsType
struct NonlinearSolveOptions <: ExternalSolverOptions

Options for the SimpleNonlinearSolve pacakge solvers

Fields

  • algorithm::SimpleNonlinearSolve algorithm = SimpleNonlinearSolve.SimpleNewtonRaphson : algorithm to use
  • additional_kwargs = (;) : any additional keyword arguments for the solver
  • atol::Float = 1e-12 : absolute convergence tolerance
  • iteration_limit::Float = 25 : maximum number of iterations
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.MinpackOptionsType
struct MinpackOptions <: ExternalSolverOptions

Options for the MINPACK's HYBRJ solver

Fields

  • algorithm::Symbol = :hybr : algorithm to use in MINPACK.jl (hybr is HYBRJ when the jacobian is provided)
  • atol::FLoat = 1e-12 : absolute convergence tolerance
  • iteration_limit::FLoat = 100 : maximum number of iterations
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.SIAMFANLEOptionsType
struct SIAMFANLEOptions <: ExternalSolverOptions

Options for the SIAMFANLEquations pacakge solvers

Fields

  • algorithm::SIAMFANLEquations algorithm = SIAMFANLEquations.nsoli : algorithm to use
  • rtol::Float = 0.0 : relative convergence tolerance
  • atol::Float = 1e-10 : absolute convergence tolerance
  • iteration_limit::Int = 1000 : maximum number of iterations
  • linear_iteration_limit::Float = 5 : maximum number of linear solve iterations (GMRES)
  • additional_kwargs = (;) : any additional keyword arguments for the solver
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.SpeedMappingOptionsType
struct SpeedMappingOptions <: ExternalSolverOptions

Options for the SpeedMapping.jl package solver

Fields

  • `orders::AbstractArray{Int} = [3, 2]
  • sig_min::Int = 0 : maybe set to 1?
  • stabilize::Bool = false : stabilizes before extrapolation
  • check_obj::Bool = false : checks for inf's and nan's and starts from previous finite point
  • atol::Float = 1e-10 : absolute convergence tolerance
  • iteration_limit::Float = 1000 : maximum number of iterations
  • time_limit::Float = Inf : time limit in seconds
  • lower::Float = nothing : box lower bounds
  • upper::Float = nothing : box upper bounds
  • buffer::Float = 0.01 : if using bounds, buffer brings x inside bounds by buffer amountd
  • Lp::Float = Inf : p value for p-norm for convergence criteria
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.FixedPointOptionsType
struct FixedPointOptions <: ExternalSolverOptions

Options for the FixedPoint.jl package solver

Fields

  • iteration_limit::Int = 1000 : maximum number of iterations
  • vel::Float = 0.9 : vel keyword argument, default is package default
  • ep::Float = 0.01 : ep keyword argument, default is package default
  • atol::Float = 1e-12 : absolute convergence tolerance
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.CSORSolverOptionsType
struct CSORSolverOptions <: InternalSolverOptions

Type containing all the options for the CSOR (controlled successive over relaxation) solver.

Note that the defaults match DFDC with the exception of the relaxation schedule, which is an experimental feature.

Fields

  • verbose::Bool = false : flag to print verbose statements
  • iteration_limit::Float = 1e2 : maximum number of iterations
  • nrf::Float = 0.4 : nominal relaxation factor
  • bt1::Float = 0.2 : backtracking factor 1
  • bt2::Float = 0.6 : backtracking factor 2
  • pf1::Float = 0.4 : press forward factor 1
  • pf2::Float = 0.5 : press forward factor 2
  • btw::Float = 0.6 : backtracking factor for wake
  • pfw::Float = 1.2 : press forward factor for wake
  • relaxation_schedule::TS = [[0.0;1e-14;1e-13;1e10]), [1.0;1.0;0.0;0.0])] : values used in spline definition for scaling the relaxation factors (second vector) after various convergence values (first vector).
  • f_circ::Float = 1e-3 : convergence tolerance for rotor circulation
  • f_dgamw::Float = 2e-4 : convergence tolerance for wake vortex strength
  • convergence_type::ConvergenceType = Relative() : dispatch for relative or absolute convergence criteria.
  • Vconv::AbstractArray{Float} = [1.0] : velocity used in relative convergence criteria (should be set to Vref).
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
source
DuctAPE.ModCSORSolverOptionsType
struct ModCSORSolverOptions <: InternalSolverOptions

Type containing all the options for the modified CSOR solver.

Fields

  • verbose::Bool = false : flag to print verbose statements
  • iteration_limit::Float = 350 : maximum number of iterations
  • relaxation_parameters::NamedTuple = (;
    • nrf::Float = 0.4 : nominal relaxation factor
    • bt1::Float = 0.2 : backtracking factor 1
    • bt2::Float = 0.6 : backtracking factor 2
    • pf1::Float = 0.4 : press forward factor 1
    • pf2::Float = 0.5 : press forward factor 2
    • btw::Float = 0.6 : backtracking factor for wake
    • pfw::Float = 1.2 : press forward factor for wake
    ) : parameters for determining relaxation level of states in each iteration.
  • converged::AbstractArray{Bool} = [false] : flag to track if convergence took place.
  • iterations::AbstractArray{Int} = [0] : iteration counter
  • residuals::AbstractArray{Int} = [0] : iteration counter
source

Preprocess

DuctAPE.setup_analysisFunction
setup_analysis(
    ducted_rotor::DuctedRotor,
    operating_point::OperatingPoint,
    options::Options=set_options();
    prepost_container_caching=nothing,
    solve_parameter_caching=nothing,
    solve_container_caching=nothing,
)

Perform pre-processing and cache setup (as needed) for propuslor analysis.

Arguments

  • ducted_rotor::DuctedRotor : DuctedRotor input object (see docstring for DuctedRotor type)
  • operating_point::OperatingPoint : OperatingPoint input object (see docstring for OperatingPoint type)
  • options::Options=set_options() : Options object (see set_options and related functions)

Keyword Arguments

  • prepost_container_caching=nothing : Output of allocate_prepost_container_cache
  • solve_parameter_caching=nothing : Output of allocate_solve_parameter_container_cache
  • solve_container_caching=nothing : Output of allocate_solve_container_cache

Returns

  • problem_dimensions::NamedTuple : Named Tuple contiaining bookkeeping information (problem dimensions)
  • prepost_containers::NamedTuple : Named Tuple containing reshaped views into the prepost cache
  • solve_parameter_cache_vector::Vector : Vector containing the relevant typed cache vector of solve parameters
  • solve_parameter_cache_dims::NamedTuple : Named Tuple containing dimensions used for reshaping the solve parameter cache
  • A_bb_LU::LinearAlgebra.LU : The LU factorization of the AIC matrix used in the panel method
  • lu_decomp_flag::Bool : flag indicating if the LU decomposition was successful
  • airfoils::Matrix{AFType} : Matrix contiaining the blade element airfoil polar objects
  • idmaps::NamedTuple : Named Tuple containing bookkeeping information (index mappings)
source

Analysis

DuctAPE.analyzeFunction
analyze(
    ducted_rotor::DuctedRotor,
    operating_point::OperatingPoint,
    reference_parameters::ReferenceParameters,
    options::Options=set_options();
    prepost_container_caching=nothing,
    solve_parameter_caching=nothing,
    solve_container_caching=nothing,
    return_inputs=false,
)

Analyze ducted_rotor, including preprocessing.

Arguments

  • ducted_rotor::DuctedRotor : DuctedRotor input object (see docstring for DuctedRotor type)
  • operating_point::OperatingPoint : OperatingPoint input object (see docstring for OperatingPoint type)
  • reference_parameters::ReferenceParameters : ReferenceParameters input object (see docstring for ReferenceParameters type)
  • options::Options=set_options() : Options object (see set_options and related functions)

Keyword Arguments

  • prepost_container_caching=nothing : Output of allocate_prepost_container_cache
  • solve_parameter_caching=nothing : Output of allocate_solve_parameter_container_cache
  • solve_container_caching=nothing : Output of allocate_solve_container_cache
  • return_inputs=false : flag as to whether or not to return the pre-processed inputs

Returns

  • outs::NamedTuple : Named Tuple of various analysis outputs (see docstring for postprocess for details), note, if linear system decomposition fails, no solve is performed and an empty vector is returned.
  • ins::NamedTuple : Named Tuple of various pre-processed inputs (e.g. panels and body linear system), will only be returned if return_inputs=true
  • convergence_flag : Flag for successful solve convergence
source
analyze(
    ducted_rotor::DuctedRotor,
    operating_point::OperatingPoint,
    reference_parameters::ReferenceParameters,
    prepost_containers,
    solve_parameter_cache_vector,
    solve_parameter_cache_dims,
    airfoils,
    A_bb_LU,
    idmaps,
    problem_dimensions,
    options::Options=set_options();
    return_inputs=false,
    solve_container_caching=nothing,
)

Analyze ducted_rotor, assuming setup_analysis has been called and the outputs thereof are being passed in here.

Arguments

  • ducted_rotor::DuctedRotor : DuctedRotor input object (see docstring for DuctedRotor type)
  • operating_point::OperatingPoint : OperatingPoint input object (see docstring for OperatingPoint type)
  • reference_parameters::ReferenceParameters : ReferenceParameters input object (see docstring for ReferenceParameters type)
  • prepost_containers::NamedTuple : An output from setup_analysis containing reshaped views into the prepost cache
  • solve_parameter_cache_vector::Vector : An output from setup_analysis containing the relevant typed cache vector of solve parameters
  • solve_parameter_cache_dims::NamedTuple : An output from setup_analysis containing dimensions used for reshaping the solve parameter cache
  • airfoils::Vector{AFType} : An output from setup_analysis contiaining the blade element airfoil polar objects
  • A_bb_LU::LinearAlgebra.LU : An output from setup_analysis that is the LU decomposition of the AIC matrix used in the panel method
  • idmaps::NamedTuple : An output from setup_analysis containing bookkeeping information (index mappings)
  • problem_dimensions::NamedTuple : An output from setup_analysis contiaining bookkeeping information (problem dimensions)
  • options::Options=set_options() : Options object

Keyword Arguments

  • solve_container_caching=nothing : Output of allocate_solve_container_cache
  • return_inputs=false : flag as to whether or not to return the pre-processed inputs

Returns

  • outs::NamedTuple : Named Tuple of various analysis outputs (see docstring for postprocess for details), note, if linear system decomposition fails, no solve is performed and an empty vector is returned.
  • ins::NamedTuple : Named Tuple of various pre-processed inputs (e.g. panels and body linear system), will only be returned if return_inputs=true
  • convergence_flag : Flag for successful solve convergence
source
analyze(
    ducted_rotor::DuctedRotor,
    operating_point::AbstractVector{OperatingPoint},
    reference_parameters::ReferenceParameters,
    options::Options=set_options();
    prepost_container_caching=nothing,
    solve_parameter_caching=nothing,
    solve_container_caching=nothing,
    return_inputs=false,
)

Analyze ducted_rotor, including preprocessing, for a set of operating points.

Arguments

  • ducted_rotor::DuctedRotor : DuctedRotor input object
  • operating_point::AbstractVector{OperatingPoint} : Vector of Operating Points at which to analyze the ducted_rotor
  • reference_parameters::ReferenceParameters : ReferenceParameters input object (see docstring for ReferenceParameters type)
  • options::Options=set_options() : Options object

Keyword Arguments

  • prepost_container_caching=nothing : Output of allocate_prepost_container_cache
  • solve_parameter_caching=nothing : Output of allocate_solve_parameter_container_cache
  • solve_container_caching=nothing : Output of allocate_solve_container_cache
  • return_inputs=false : flag as to whether or not to return the pre-processed inputs

Returns

  • outs::Vector{NamedTuple} : Vector of named tuples of various analysis outputs (see docstring for postprocess for details), note, if linear system decomposition fails, no solve is performed and an empty vector is returned.
  • ins::NamedTuple : Named Tuple of various pre-processed inputs (e.g. panels and body linear system), will only be returned if return_inputs=true
  • convergence_flag : Flag for successful solve convergence
source
analyze(
    ducted_rotor::DuctedRotor,
    operating_point::Vector{OperatingPoint},
    reference_parameters::ReferenceParameters,
    prepost_containers,
    solve_parameter_cache_vector,
    solve_parameter_cache_dims,
    airfoils,
    A_bb_LU,
    idmaps,
    problem_dimensions,
    options::Options=set_options();
    return_inputs=false,
    solve_container_caching=nothing,
)

Analyze ducted_rotor, assuming setup_analysis has been called and the inputs are being passed in here.

Arguments

  • ducted_rotor::DuctedRotor : DuctedRotor input object
  • operating_point::AbstractVector{OperatingPoint} : Vector of Operating Points at which to analyze the ducted_rotor
  • reference_parameters::ReferenceParameters : ReferenceParameters input object (see docstring for ReferenceParameters type)
  • prepost_containers::NamedTuple : An output from setup_analysis containing reshaped views into the prepost cache
  • solve_parameter_cache_vector::Vector : An output from setup_analysis containing the relevant typed cache vector of solve parameters
  • solve_parameter_cache_dims::NamedTuple : An output from setup_analysis containing dimensions used for reshaping the solve parameter cache
  • airfoils::Vector{AFType} : An output from setup_analysis contiaining the blade element airfoil polar objects
  • A_bb_LU::LinearAlgebra.LU : An output from setup_analysis that is the LU decomposition of the AIC matrix used in the panel method
  • idmaps::NamedTuple : An output from setup_analysis containing bookkeeping information (index mappings)
  • problem_dimensions::NamedTuple : An output from setup_analysis contiaining bookkeeping information (problem dimensions)
  • options::Options=set_options() : Options object

Keyword Arguments

  • solve_container_caching=nothing : Output of allocate_solve_container_cache
  • return_inputs=false : flag as to whether or not to return the pre-processed inputs

Returns

  • outs::Vector{NamedTuple} : Named Tuple of various analysis outputs (see docstring for postprocess for details), note, if linear system decomposition fails, no solve is performed and an empty vector is returned.
  • ins::NamedTuple : Named Tuple of various pre-processed inputs (e.g. panels and body linear system), will only be returned if return_inputs=true. Note that some inputs will be overwritten (e.g. the linear system RHS components related to the freestream) and only those associated with the final operating point will be returned.
  • convergence_flag : Flag for successful solve convergence
source

Postprocess

DuctAPE.HeadsBoundaryLayerOptionsType
struct HeadsBoundaryLayerOptions

Fields:

  • model_drag::Tb=false : flag to turn on viscous drag approximation
  • n_steps::Int = Int(5e2) : number of steps to use in boundary layer integration
  • first_step_size::Float = 1e-6 : size of first step in boundary layer integration
  • offset::Float = 1e-3 : size of offset for (where to initialize) boundary layer integration
  • solver_type::Type = DiffEq() : type of ODE solver (RK() or DiffEq())
  • ode::Function = RadauIIA5 : solver to use for boundary layer integration (RadauIIA5, RK4, or RK2 available)
  • separation_criteria::Float=3.0 : value of H12 after which separation should happen.
  • separation_allowance_upper::Int=10 : upper side allowance for how many steps ahead of the trailing edge we'll allow separation without penalty
  • separation_allowance_lower::Int=10 : lower side allowance for how many steps ahead of the trailing edge we'll allow separation without penalty
  • separation_penalty_upper::Float=0.2 : upper side maximum penalty value for separation (at leading edge)
  • separation_penalty_lower::Float=0.2 : lower side maximum penalty value for separation (at leading edge)
source
DuctAPE.GreensBoundaryLayerOptionsType
struct GreensBoundaryLayerOptions

NOTE: Green's method is mostly implemented, but there are several bugs still, especially when using Imperial units. Known Bugs:

  • Imperial units overestimate momentum thickness. Likely a unit conversion bug.
  • In some cases of non-separation, the momentum thickens or shape parameter becomes exceedingly large, vastly overestimating the drag coefficient. Likely the product of one or more of the adjustments to try and make the method more robust.

Fields:

  • model_drag::Tb=true : flag to turn off viscous drag approximation
  • lambda::Bool = true : flag to add secondary influences into boundary layer residuals
  • longitudinal_curvature::Bool = true : if lambda=true, flag to add longitudinal curvature influence into boundary layer residuals
  • lateral_strain::Bool = true : if lambda=true, flag to add lateral strain influence into boundary layer residuals
  • dilation::Bool = true : if lambda=true, flag to add dilation influence into boundary layer residuals
  • n_steps::Int = Int(2e2) : number of steps to use in boundary layer integration
  • first_step_size::Float = 1e-3 : size of first step in boundary layer integration
  • offset::Float = 1e-2 : size of offset for (where to initialize) boundary layer integration
  • solver_type::Type = DiffEq() : type of ODE solver (RK() or DiffEq())
  • ode::Function = RadauIIA5 : solver to use for boundary layer integration (RadauIIA5, RK4, or RK2 available)
  • separation_allowance_upper::Int=3 : upper side allowance for how many steps ahead of the trailing edge we'll allow separation without penalty
  • separation_allowance_lower::Int=3 : lower side allowance for how many steps ahead of the trailing edge we'll allow separation without penalty
  • separation_penalty_upper::Float=0.2 : upper side maximum penalty value for separation (at leading edge)
  • separation_penalty_lower::Float=0.2 : lower side maximum penalty value for separation (at leading edge)
source

Miscellaneous

Airfoil/Geometry Manipulation

NACA 6-Series Cascade Geometry Generation