Public API
Input Types
DuctAPE.PanelingConstants
— TypePanelingConstants(
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.
DuctAPE.Rotor
— TypeRotor(
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 selectedrotorzloc
.Rtip::AbstractVector{Float}
: Dimensional tip radius of rotor. Is used to determine the radial position of the duct if theautoshiftduct
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 withC4Blade.DFDCairfoil
types.fliplift::AbstractVector{Bool}
: Flag to indicate if the airfoil lift values should be flipped or not.
DuctAPE.DuctedRotor
— TypeDuctedRotor(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 theautoshiftduct
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.
DuctAPE.OperatingPoint
— TypeOperatingPoint(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 numberrhoinf::AbstractVector{Float}
: Freestream densitymuinf::AbstractVector{Float}
: Freestream viscosityasound::AbstractVector{Float}
: Freestream speed of soundPtot::AbstractVector{Float}
: Freestream total pressureTtot::AbstractVector{Float}
: Freestream total temperatureOmega::AbstractVector{Float}
: Rotor rototation rate(s)
DuctAPE.ReferenceParameters
— TypeReferenceParameters(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.
Preallocations
DuctAPE.allocate_prepost_container_cache
— Functionallocate_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 cacheprepost_container_cache_dims::NamedTuple
: a named tuple containing the dimensions used for reshaping the cache when needed.
DuctAPE.allocate_solve_parameter_cache
— Functionallocate_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 dispatchpaneling_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 cachesolve_parameter_cache_dims::NamedTuple
: a named tuple containing the dimensions used for reshaping the cache when needed.
DuctAPE.allocate_solve_container_cache
— Functionallocate_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 cachesolve_container_cache_dims::NamedTuple
: a named tuple containing the dimensions used for reshaping the cache when needed.
Options
General Options
DuctAPE.Options
— Typestruct Options
Type containing (nearly) all the available user options.
Fields
General Options
verbose::Bool = false
: flag to print verbose statementssilence_warnings::Bool = true
: flag to silence warningsmultipoint_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 bodiesautoshiftduct::Bool = true
: flag as to whether duct geometry should be shifted based on rotor tip locationlu_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 axistegaptol::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 objectwrite_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 writecheckoutfileexists::Bool = false
: Flag for whether to check if file exists before overwritingoutput_tuple_name::AbstractArray{String} = ["outs"]
: variable name for named tuple written to out file
Solving Options
grid_solver_options::GridSolverOptionsType = GridSolverOptions()
: elliptic grid solver optionssolver_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.
DuctAPE.set_options
— Functionset_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.
Integration Options
DuctAPE.IntegrationOptions
— Typestruct 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.
DuctAPE.GaussLegendre
— Typestruct GaussLegendre <: IntegrationMethod
Options for Gauss-Legendre integration method
Fields
sample_points::Int
: Sample Pointsweights::Int
: Gauss weights
DuctAPE.GaussKronrod
— Typestruct GaussKronrod <: IntegrationMethod
Options for Gauss-Kronrod integration method
Fields
order::Int = 7
: order of Legendre polynomial to use on each intervalmaxevales::Int = 10^7
: maximum number of evaluations in the adaptive methodatol::Float = 0.0
: absolute error tolerance. (note, if zero, QuadGK uses sqrt(eps()) relative tolerance).
DuctAPE.Romberg
— Typestruct 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.
Solver Options
Elliptic Grid Solve
DuctAPE.SLORGridSolverOptions
— Typestruct SLORGridSolverOptions <: GridSolverOptionsType
Options for SLOR (successive line over relaxation) elliptic grid solver.
Fields
iteration_limit::Int = 100
: maximum number of iterationsatol::Float = 1e-9
: absolute convergence toleranceconverged::AbstractArray{Bool}
= [false]iterations::AbstractArray{Int} = [0]
: iteration counter
DuctAPE.GridSolverOptions
— Typestruct GridSolverOptions <: GridSolverOptionsType
Options for Newton elliptic grid solver.
Fields
iteration_limit::Int = 20
: maximum number of iterationsatol::Float = 3e-10
: absolute convergence tolerancealgorithm::Symbol = :newton
: algorithm to use in NLsolve.jlautodiff::Symbol = :forward
: differentiation method to use in NLsolve.jlprecondition = false
: flag to precondition with SLORprecondition_max_iterations = 3
: number of precondition iterationsconverged::AbstractArray{Bool}
= [false]iterations::AbstractArray{Int} = [0]
: iteration counterresidual_value::AbstractArray{Int} = [0]
: residual value
Aerodynamics Solve
DuctAPE.ChainSolverOptions
— Typestruct 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
DuctAPE.CompositeSolverOptions
— Typestruct 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
DuctAPE.NLsolveOptions
— Typestruct NLsolveOptions <: ExternalSolverOptions
Options for the NLsolve pacakge solvers
Fields
algorithm::Symbol = :anderson
: algorithm to useadditional_kwargs = (;)
: any additional keyword arguments for the solveratol::Float = 1e-12
: absolute convergence toleranceiteration_limit::Int = 25
: maximum number of iterationslinesearch_method::LineSearches method = LineSearches.MoreThuente
: line search method to uselinesearch_kwargs = (;)
: any additional lineseach keyword argumentsconverged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counter
DuctAPE.NonlinearSolveOptions
— Typestruct NonlinearSolveOptions <: ExternalSolverOptions
Options for the SimpleNonlinearSolve pacakge solvers
Fields
algorithm::SimpleNonlinearSolve algorithm = SimpleNonlinearSolve.SimpleNewtonRaphson
: algorithm to useadditional_kwargs = (;)
: any additional keyword arguments for the solveratol::Float = 1e-12
: absolute convergence toleranceiteration_limit::Float = 25
: maximum number of iterationsconverged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counter
DuctAPE.MinpackOptions
— Typestruct 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 toleranceiteration_limit::FLoat = 100
: maximum number of iterationsconverged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counter
DuctAPE.SIAMFANLEOptions
— Typestruct SIAMFANLEOptions <: ExternalSolverOptions
Options for the SIAMFANLEquations pacakge solvers
Fields
algorithm::SIAMFANLEquations algorithm = SIAMFANLEquations.nsoli
: algorithm to usertol::Float = 0.0
: relative convergence toleranceatol::Float = 1e-10
: absolute convergence toleranceiteration_limit::Int = 1000
: maximum number of iterationslinear_iteration_limit::Float = 5
: maximum number of linear solve iterations (GMRES)additional_kwargs = (;)
: any additional keyword arguments for the solverconverged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counter
DuctAPE.SpeedMappingOptions
— Typestruct 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 extrapolationcheck_obj::Bool = false
: checks for inf's and nan's and starts from previous finite pointatol::Float = 1e-10
: absolute convergence toleranceiteration_limit::Float = 1000
: maximum number of iterationstime_limit::Float = Inf
: time limit in secondslower::Float = nothing
: box lower boundsupper::Float = nothing
: box upper boundsbuffer::Float = 0.01
: if using bounds, buffer brings x inside bounds by buffer amountdLp::Float = Inf
: p value for p-norm for convergence criteriaconverged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counter
DuctAPE.FixedPointOptions
— Typestruct FixedPointOptions <: ExternalSolverOptions
Options for the FixedPoint.jl package solver
Fields
iteration_limit::Int = 1000
: maximum number of iterationsvel::Float = 0.9
: vel keyword argument, default is package defaultep::Float = 0.01
: ep keyword argument, default is package defaultatol::Float = 1e-12
: absolute convergence toleranceconverged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counter
DuctAPE.CSORSolverOptions
— Typestruct 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 statementsiteration_limit::Float = 1e2
: maximum number of iterationsnrf::Float = 0.4
: nominal relaxation factorbt1::Float = 0.2
: backtracking factor 1bt2::Float = 0.6
: backtracking factor 2pf1::Float = 0.4
: press forward factor 1pf2::Float = 0.5
: press forward factor 2btw::Float = 0.6
: backtracking factor for wakepfw::Float = 1.2
: press forward factor for wakerelaxation_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 circulationf_dgamw::Float = 2e-4
: convergence tolerance for wake vortex strengthconvergence_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
DuctAPE.ModCSORSolverOptions
— Typestruct ModCSORSolverOptions <: InternalSolverOptions
Type containing all the options for the modified CSOR solver.
Fields
verbose::Bool = false
: flag to print verbose statementsiteration_limit::Float = 350
: maximum number of iterationsrelaxation_parameters::NamedTuple
= (;nrf::Float = 0.4
: nominal relaxation factorbt1::Float = 0.2
: backtracking factor 1bt2::Float = 0.6
: backtracking factor 2pf1::Float = 0.4
: press forward factor 1pf2::Float = 0.5
: press forward factor 2btw::Float = 0.6
: backtracking factor for wakepfw::Float = 1.2
: press forward factor for wake
converged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counterresiduals::AbstractArray{Int} = [0]
: iteration counter
Preprocess
DuctAPE.setup_analysis
— Functionsetup_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 forDuctedRotor
type)operating_point::OperatingPoint
: OperatingPoint input object (see docstring forOperatingPoint
type)options::Options=set_options()
: Options object (seeset_options
and related functions)
Keyword Arguments
prepost_container_caching=nothing
: Output ofallocate_prepost_container_cache
solve_parameter_caching=nothing
: Output ofallocate_solve_parameter_container_cache
solve_container_caching=nothing
: Output ofallocate_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 cachesolve_parameter_cache_vector::Vector
: Vector containing the relevant typed cache vector of solve parameterssolve_parameter_cache_dims::NamedTuple
: Named Tuple containing dimensions used for reshaping the solve parameter cacheA_bb_LU::LinearAlgebra.LU
: The LU factorization of the AIC matrix used in the panel methodlu_decomp_flag::Bool
: flag indicating if the LU decomposition was successfulairfoils::Matrix{AFType}
: Matrix contiaining the blade element airfoil polar objectsidmaps::NamedTuple
: Named Tuple containing bookkeeping information (index mappings)
Analysis
DuctAPE.analyze
— Functionanalyze(
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 forDuctedRotor
type)operating_point::OperatingPoint
: OperatingPoint input object (see docstring forOperatingPoint
type)reference_parameters::ReferenceParameters
: ReferenceParameters input object (see docstring forReferenceParameters
type)options::Options=set_options()
: Options object (seeset_options
and related functions)
Keyword Arguments
prepost_container_caching=nothing
: Output ofallocate_prepost_container_cache
solve_parameter_caching=nothing
: Output ofallocate_solve_parameter_container_cache
solve_container_caching=nothing
: Output ofallocate_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 ifreturn_inputs=true
convergence_flag
: Flag for successful solve convergence
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 forDuctedRotor
type)operating_point::OperatingPoint
: OperatingPoint input object (see docstring forOperatingPoint
type)reference_parameters::ReferenceParameters
: ReferenceParameters input object (see docstring forReferenceParameters
type)prepost_containers::NamedTuple
: An output fromsetup_analysis
containing reshaped views into the prepost cachesolve_parameter_cache_vector::Vector
: An output fromsetup_analysis
containing the relevant typed cache vector of solve parameterssolve_parameter_cache_dims::NamedTuple
: An output fromsetup_analysis
containing dimensions used for reshaping the solve parameter cacheairfoils::Vector{AFType}
: An output fromsetup_analysis
contiaining the blade element airfoil polar objectsA_bb_LU::LinearAlgebra.LU
: An output fromsetup_analysis
that is the LU decomposition of the AIC matrix used in the panel methodidmaps::NamedTuple
: An output fromsetup_analysis
containing bookkeeping information (index mappings)problem_dimensions::NamedTuple
: An output fromsetup_analysis
contiaining bookkeeping information (problem dimensions)options::Options=set_options()
: Options object
Keyword Arguments
solve_container_caching=nothing
: Output ofallocate_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 ifreturn_inputs=true
convergence_flag
: Flag for successful solve convergence
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 objectoperating_point::AbstractVector{OperatingPoint}
: Vector of Operating Points at which to analyze the ducted_rotorreference_parameters::ReferenceParameters
: ReferenceParameters input object (see docstring forReferenceParameters
type)options::Options=set_options()
: Options object
Keyword Arguments
prepost_container_caching=nothing
: Output ofallocate_prepost_container_cache
solve_parameter_caching=nothing
: Output ofallocate_solve_parameter_container_cache
solve_container_caching=nothing
: Output ofallocate_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 ifreturn_inputs=true
convergence_flag
: Flag for successful solve convergence
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 objectoperating_point::AbstractVector{OperatingPoint}
: Vector of Operating Points at which to analyze the ducted_rotorreference_parameters::ReferenceParameters
: ReferenceParameters input object (see docstring forReferenceParameters
type)prepost_containers::NamedTuple
: An output fromsetup_analysis
containing reshaped views into the prepost cachesolve_parameter_cache_vector::Vector
: An output fromsetup_analysis
containing the relevant typed cache vector of solve parameterssolve_parameter_cache_dims::NamedTuple
: An output fromsetup_analysis
containing dimensions used for reshaping the solve parameter cacheairfoils::Vector{AFType}
: An output fromsetup_analysis
contiaining the blade element airfoil polar objectsA_bb_LU::LinearAlgebra.LU
: An output fromsetup_analysis
that is the LU decomposition of the AIC matrix used in the panel methodidmaps::NamedTuple
: An output fromsetup_analysis
containing bookkeeping information (index mappings)problem_dimensions::NamedTuple
: An output fromsetup_analysis
contiaining bookkeeping information (problem dimensions)options::Options=set_options()
: Options object
Keyword Arguments
solve_container_caching=nothing
: Output ofallocate_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 ifreturn_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
Postprocess
DuctAPE.BoundaryLayerOptions
— Typeabstract type BoundaryLayerOptions
Used in boundary layer method dispatch
DuctAPE.HeadsBoundaryLayerOptions
— Typestruct HeadsBoundaryLayerOptions
Fields:
model_drag::Tb=false
: flag to turn on viscous drag approximationn_steps::Int = Int(5e2)
: number of steps to use in boundary layer integrationfirst_step_size::Float = 1e-6
: size of first step in boundary layer integrationoffset::Float = 1e-3
: size of offset for (where to initialize) boundary layer integrationsolver_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 penaltyseparation_allowance_lower::Int=10
: lower side allowance for how many steps ahead of the trailing edge we'll allow separation without penaltyseparation_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)
DuctAPE.GreensBoundaryLayerOptions
— Typestruct 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 approximationlambda::Bool = true
: flag to add secondary influences into boundary layer residualslongitudinal_curvature::Bool = true
: iflambda
=true, flag to add longitudinal curvature influence into boundary layer residualslateral_strain::Bool = true
: iflambda
=true, flag to add lateral strain influence into boundary layer residualsdilation::Bool = true
: iflambda
=true, flag to add dilation influence into boundary layer residualsn_steps::Int = Int(2e2)
: number of steps to use in boundary layer integrationfirst_step_size::Float = 1e-3
: size of first step in boundary layer integrationoffset::Float = 1e-2
: size of offset for (where to initialize) boundary layer integrationsolver_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 penaltyseparation_allowance_lower::Int=3
: lower side allowance for how many steps ahead of the trailing edge we'll allow separation without penaltyseparation_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)