Options
There are quite a few options to choose from in DuctAPE. As DuctAPE was developed, various options were added as different approaches were considered. The current defaults were selected as the best set of options for the studies in the publications metioned in DuctAPE's README. However, most of the options implemented during development have been maintained in case of future need. The Options
object contains various general options as well as several more detailed option objects. To make the process of setting options easier, the set_options
method can be used to make only the desired changes from the defaults without having to define everything in the Options
object. Note that the Options
object is implemented using the @kwdef
macro, so the set_options
function doesn't really do anything in the case of a single operating point, but the other dispatch of set_options
is especially helpful in initializing several of the sub-option objects to defaults of the correct size. Note that there are several fields in the options that are used for bookkeeping, especially in the case of multiple operating points.
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 warningshard_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.
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 location
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
Bookkeeping Options
multipoint_index::Int = [1]
: holds current index of multi-point solver (no need for user to change this usually)lu_decomp_flag::Bool = false
: flag indicating if panel method LHS matrix factorization was successful
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.
The major sub-categories of options include general options, pre-processing options, solver options for both determining the wake sheet positions as well as the overall aerodyanmics solve, post-processing options, and bookkeeping options.
Bookkeeping Options
These are options that can be changed by the user for development/debugging purposes, but at this point, it would be wise in general usage to not change them. In future revisions, these will likely no longer be accessible to the user.
General Options
The verbose and silence warnings options are simply about what get's printed as the analysis runs. Warnings are printed when some sort of automated adjustment is made to the inputs in order to ensure they conform to the format required. The verbose option at this level is for verbose statements that are not within any solvers. Solver verbosity is constrolled in the individual solver options.
Occasionally, something in the preprocessing will fail, likely the LU decomposition of the linear system defining the bodies' panel system. If such a failure occurs, DuctAPE cannot continue to the main solve and will exit. The hard_fail
option dictates what the exit behavior is. If true, DuctAPE will just return nothing
immediately, which is quicker for turn-around on single runs. If false, DuctAPE will attempt to return an output object of the correct size and type, which is convenient for some optimization frameworks for which you'll want some output to be available even if passing a failure flag for the specific analysis.
Preprocess Options
Geometry Interpolation and Generation Options
The autoshiftduct
option may be convenient depending on how the duct coordinates are being input. It allows the user to input the duct coordinates at an arbitrary radial location, for example if a standard airfoil is used with leading edge at (0,0). In the preprocessing, the duct geometry will be shifted to the radial location at which the rotor tip is coincident with the duct surface at the axial location at which the first rotor is situated. If you are already inputting the duct geometry at the correct position, this option may be turned off, but it usually doesn't hurt to be left on.
Paneling Options
These, in general, do not need to be touched by users; thus we do not include them in the PanelingConstants
, but we do make them available.
Integration Options
DuctAPE uses numerical integration to determine the influence of axisymmetric vortex and source panels on the restof the system. There are several options for numerical integration, and the methods can be mixed and matched with their own specific options for the nominal (panel on other panels) and singular (panel on itself) cases. The Gauss-Legendre method is useful in optimization cases to avoid noise in the integration error. The Gauss-Kronrod method is implemented via QuadGK.jl and is more accurate, especially in for the singular cases, but is less useful for optimization purposes due to noise in the error from the adaptive nature. Similarly, the Romberg method similar to that implemented in DFDC is adaptive and can be fast, but is also not usually the best choice. Thus we default to Gauss-Legendre methods.
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.
Example
using DuctAPE
# set nominal options using a GaussLegendre object (which is an InterationMethod type)
# note that a convenience method is used here that takes in the number of points and
#calculates the appropriate sample locations and weights.
nominal_integration_method = DuctAPE.GaussLegendre(10)
# set singular options using a GaussKronrod object (which is an InterationMethod type)
# note that like most option structs, these are defined using @kwdef allowing the fields
#to be treated as keyword arguments.
# also note that we haven't changed the evaluation limit (default 10^7)
singular_integration_method = DuctAPE.GaussKronrod(; order=7, atol=2e-16)
# put the quadrature options together
integration_options = DuctAPE.IntegrationOptions(;
nominal=nominal_integration_method, singular=singular_integration_method
)
# example of calling the set_options function
options = DuctAPE.set_options(; integration_options=integration_options)
Options{Bool, DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}, Vector{Bool}, Float64, Vector{Int64}, Vector{String}, Vector{String}, DuctAPE.var"#56#61", IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussKronrod{Float64, Int64}}, ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}, GridSolverOptions{Bool, Float64, Int64, Symbol}}(false, true, [1], DuctAPE.var"#56#61"(), true, false, 0.05, 1.0e-15, 2.220446049250313e-15, IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussKronrod{Float64, Int64}}(GaussLegendre{Vector{Float64}, Vector{Float64}}([0.013046735741414128, 0.06746831665550773, 0.16029521585048778, 0.2833023029353764, 0.4255628305091844, 0.5744371694908156, 0.7166976970646236, 0.8397047841495122, 0.9325316833444923, 0.9869532642585859], [0.033335672154344104, 0.07472567457529028, 0.10954318125799103, 0.13463335965499826, 0.14776211235737646, 0.14776211235737646, 0.13463335965499826, 0.10954318125799103, 0.07472567457529028, 0.033335672154344104]), GaussKronrod{Float64, Int64}(7, 10000000, 2.0e-16)), DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}(false, true, true, false, 5000, 1.0e-6, nothing, nothing, 0.001, DuctAPE.RK(), DuctAPE.RK2, 3.0, 10, 10, 0.2, 0.2, false, 0.0, 0.0001, 0.0, false), Bool[0], ["outputs.jl"], false, ["outs"], GridSolverOptions{Bool, Float64, Int64, Symbol}(30, 3.0e-10, :newton, :forward, false, 3, Bool[0], [0], [0.0]), ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}(false, 500, (nrf = 0.4, bt1 = 0.2, bt2 = 0.6, pf1 = 0.4, pf2 = 0.5, btw = 0.6, pfw = 1.2), 5.0e-10, Bool[0], [0], [-1.0]), true)
Solver Options
Elliptic Grid Solver Options
The wake geometry is obtained by solving for approximate streamlines on an elliptic grid with the bodies as boundaries. There are two methods available, a successive line over relaxiation (SLOR) method that can be used in isolation or as a preconditioner to a Newton solve method from the NLsolve.jl package. The default option is to run a few iterations of SLOR to precondition and smooth out the initial grid, and then finish up with the Newton solve, usually within 3-5 iterations. Note that the SLOR method is not implemented using ImplicitAD.jl since there isn't a clean residual definition separate from the solve method. The Newton solve is implemented using ImplicitAD.jl which helps speed up automatic differentiation in an optimization setting, but it still benefits from a few iterations of SLOR beforehand.
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
Example
using DuctAPE
# define wake grid solver settings
wake_solve_options = DuctAPE.GridSolverOptions(; atol=1e-10)
# set all options
options = DuctAPE.set_options(; grid_solver_options=wake_solve_options)
Options{Bool, DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}, Vector{Bool}, Float64, Vector{Int64}, Vector{String}, Vector{String}, DuctAPE.var"#56#61", IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}, ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}, GridSolverOptions{Bool, Float64, Int64, Symbol}}(false, true, [1], DuctAPE.var"#56#61"(), true, false, 0.05, 1.0e-15, 2.220446049250313e-15, IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}(GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838]), GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838])), DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}(false, true, true, false, 5000, 1.0e-6, nothing, nothing, 0.001, DuctAPE.RK(), DuctAPE.RK2, 3.0, 10, 10, 0.2, 0.2, false, 0.0, 0.0001, 0.0, false), Bool[0], ["outputs.jl"], false, ["outs"], GridSolverOptions{Bool, Float64, Int64, Symbol}(30, 1.0e-10, :newton, :forward, false, 3, Bool[0], [0], [0.0]), ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}(false, 500, (nrf = 0.4, bt1 = 0.2, bt2 = 0.6, pf1 = 0.4, pf2 = 0.5, btw = 0.6, pfw = 1.2), 5.0e-10, Bool[0], [0], [-1.0]), true)
Aerodynamics Solver Options
Quite a few solve methods were explored in the development of DuctAPE which can be separated into three broad categories: Fixed-point iteration solvers, quasi-Newton solvers, and Newton solvers. In general, the fixed-point solvers have been faster and more robust than other methods, but all the methods have been kept in case they are desired for future development.
There are two methods that are implemented directly in DuctAPE: the CSOR and ModCSOR methods which are the controlled successive over relaxation fixed-point approach taken in DFDC and a modified version compatible with ImplicitAD.jl that is the current default and is currently the best (fastest/most robust) for optimization.
DuctAPE.CSORSolverOptions
— Typestruct CSORSolverOptions <: InternalSolverOptions
Type containing all the options for the CSOR (controlled successive over relaxation) solver.
Note that the defaults match DFDC settings.
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 wakef_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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
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{Float} = [-1.0]
: final residual values
The other methods are implemented via external dependencies and some do better than others. In our experience, the other strictly fixed-point methods work relatively well, but are middle of the road.
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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
The quasi-Newton methods are hit or miss, with Minpack doing well enough to make it into some of the compound solver options discussed below, but we have had very little success wth SIAMFANLEquations up to this point.
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-10
: 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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
NonlinearSolve is generally faster than NLsolve if the problem is large enough, but we find it to be significantly less robust. NLsolve's Anderson method is perhaps the best external method in terms of speed and robustness, but is just barely edged out by the CSOR methods.
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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
Finally, there are several compound solve methods implemented, the first chaining solvers together. The solver chain can be defined as the user wishes, but the defaults start with the fixed-point Anderson method, move to the Minpack quasi-Newton method if the fixed-point method doesn't converge in the given number of iterations, and then finishes with a full Newton method if the quasi-newton method doesn't converge. The other compound solver combines solvers in a composite manner, typically starting with a few iterations of a full Newton method to get the solver going in the right direction and then finishing with a fixed-point method. The ChainSolverOptions
was at one point the default method, but once the modified CSOR method was developed, the compound solvers weren't used much. Note that due to the way the solvers are implemented and dispatched, it is currently not possible to mix and match the CSOR methods with any of the external package methods.
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-10), MinpackOptions(; atol=1e-10), 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 counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
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-10), ]
: Vector of solver options to use.converged::AbstractArray{Bool} = [false]
: flag to track if convergence took place.iterations::AbstractArray{Int} = [0]
: iteration counterresiduals::AbstractArray{Float} = [-1.0]
: final residual values
Example
using DuctAPE
using LineSearches
# Define settings for NLsolve's newton method
aero_solver_options = DuctAPE.NLsolveOptions(;
algorithm=:newton,
atol=1e-10,
iteration_limit=30,
linesearch_method=LineSearches.BackTracking, #don't include parentheses on method handle
linesearch_kwargs=(; order=3, maxstep=1e6),
)
# set all the options
DuctAPE.set_options(; solver_options=aero_solver_options)
Options{Bool, DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}, Vector{Bool}, Float64, Vector{Int64}, Vector{String}, Vector{String}, DuctAPE.var"#56#61", IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}, NLsolveOptions{Bool, Float64, Int64, UnionAll, @NamedTuple{order::Int64, maxstep::Float64}, Symbol}, GridSolverOptions{Bool, Float64, Int64, Symbol}}(false, true, [1], DuctAPE.var"#56#61"(), true, false, 0.05, 1.0e-15, 2.220446049250313e-15, IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}(GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838]), GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838])), DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}(false, true, true, false, 5000, 1.0e-6, nothing, nothing, 0.001, DuctAPE.RK(), DuctAPE.RK2, 3.0, 10, 10, 0.2, 0.2, false, 0.0, 0.0001, 0.0, false), Bool[0], ["outputs.jl"], false, ["outs"], GridSolverOptions{Bool, Float64, Int64, Symbol}(30, 3.0e-10, :newton, :forward, false, 3, Bool[0], [0], [0.0]), NLsolveOptions{Bool, Float64, Int64, UnionAll, @NamedTuple{order::Int64, maxstep::Float64}, Symbol}(:newton, 1.0e-10, 30, LineSearches.BackTracking, (order = 3, maxstep = 1.0e6), Bool[0], [0], [-1.0]), true)
The iterations
field (not to be confused with the iterations_limit
field) in the solver options should generally not be changed. They automatically save (in-place) the number of iterations the solver performs and can be accessed after the analysis is run.
Postprocess Options
Most of the postprocess options have to do with writing the outputs to files, with the default behavior being to not write anything. These options can be useful for debugging purposes or saving outputs, though it is usually more efficient to manually save a select few outputs rather than all the outputs. The one option that is slightly more involved is the boundary layer options. Currently, only Head's method is fully implemented, but there has also been some development started on Green's method. The boundary layer options include choices regarding type of solver, with a simple 2nd-order Runge-Kutta method being the default (and appears to be best for optimization). There is also a 4th-order Runge-Kutta method implemented as well as the RadauIIA5 method from DifferentialEquations.jl which may be more accurate for single runs. Note that the default setting is to not run the boundary layer method. The model_drag
option in the boundary layer options needs to be set to true if it is desired to include the drag model.
DuctAPE.HeadsBoundaryLayerOptions
— Typestruct HeadsBoundaryLayerOptions
Fields:
model_drag::Bool=false
: flag to turn on viscous drag approximationtermiante::Bool=true
: flag to terminate solver when separation criteria is met.return_last_max_shape_factor::Bool=true
: return the last maximum shape factor to avoid false drops near the trailing edge.cutoff_Hsep::Bool=false
: Cutoff returned Hsep vector at separation criterian_steps::Int = Int(5e3)
: number of steps to use in boundary layer integrationfirst_step_size::Float = 1e-6
: size of first step in boundary layer integrationupper_step_size::Float=nothing
: uses fixed step size rather than total number of steps.lower_step_size::Float=nothing
: uses fixed step size rather than total number of steps.offset::Float = 1e-3
: size of offset for (where to initialize) boundary layer integrationsolver_type::AbstractODESolverType=RK()
: type of ODE solver (RK() or DiffEq())ode::Function=RK2
: 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)apply_separation_penalty_to_rotor::Bool=false
: flag to apply separation penalty to rotor performance.dy_eps::Float=0.0
: temporary development parameter.H1_eps::Float=1e-4
: temporary development parameter.H_eps::Float=0.0
: temporary development parameter.verbose::Bool=false
: flag to print verbose statements each iterations (beware; it's a lot)
Example
using DuctAPE
# Define Boundary Layer Settings
boundary_layer_options = DuctAPE.HeadsBoundaryLayerOptions(;
model_drag=true,
separation_penalty_upper=0.1,
separation_penalty_lower=0.1,
separation_allowance_upper=3,
separation_allowance_lower=25,
)
# set all the options
DuctAPE.set_options(; boundary_layer_options=boundary_layer_options)
Options{Bool, DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}, Vector{Bool}, Float64, Vector{Int64}, Vector{String}, Vector{String}, DuctAPE.var"#56#61", IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}, ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}, GridSolverOptions{Bool, Float64, Int64, Symbol}}(false, true, [1], DuctAPE.var"#56#61"(), true, false, 0.05, 1.0e-15, 2.220446049250313e-15, IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}(GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838]), GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838])), DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}(true, true, true, false, 5000, 1.0e-6, nothing, nothing, 0.001, DuctAPE.RK(), DuctAPE.RK2, 3.0, 3, 25, 0.1, 0.1, false, 0.0, 0.0001, 0.0, false), Bool[0], ["outputs.jl"], false, ["outs"], GridSolverOptions{Bool, Float64, Int64, Symbol}(30, 3.0e-10, :newton, :forward, false, 3, Bool[0], [0], [0.0]), ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}(false, 500, (nrf = 0.4, bt1 = 0.2, bt2 = 0.6, pf1 = 0.4, pf2 = 0.5, btw = 0.6, pfw = 1.2), 5.0e-10, Bool[0], [0], [-1.0]), true)
Advanced Options for Multi-point analyses
For using advanced options in multi-point analyses, there are various changes that need to be made to avoid run-time errors. Here is an example for setting options with the CSOR solver.
using DuctAPE
# number of operating points to analyze
nop = 3
options = DuctAPE.set_options(;
solver_options=DuctAPE.ModCSORSolverOptions(;
converged=fill(false, (1, nop)), # need a convergence flag for each operating point
iterations=zeros(Int, (1, nop)), # need a iteration count for each operating point
),
write_outputs=fill(false, nop), # we need to know which of the operating point outputs to write
outfile=fill("", nop), # we need to include names, even if they won't be used.
output_tuple_name=fill("outs", nop), # we need to include names, even if they won't be used.
)
Options{Bool, DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}, Vector{Bool}, Float64, Vector{Int64}, Vector{String}, Vector{String}, DuctAPE.var"#56#61", IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}, ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}, GridSolverOptions{Bool, Float64, Int64, Symbol}}(false, true, [1], DuctAPE.var"#56#61"(), true, false, 0.05, 1.0e-15, 2.220446049250313e-15, IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}(GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838]), GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838])), DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}(false, true, true, false, 5000, 1.0e-6, nothing, nothing, 0.001, DuctAPE.RK(), DuctAPE.RK2, 3.0, 10, 10, 0.2, 0.2, false, 0.0, 0.0001, 0.0, false), Bool[0, 0, 0], ["", "", ""], false, ["outs", "outs", "outs"], GridSolverOptions{Bool, Float64, Int64, Symbol}(30, 3.0e-10, :newton, :forward, false, 3, Bool[0], [0], [0.0]), ModCSORSolverOptions{Bool, Float64, Int64, @NamedTuple{nrf::Float64, bt1::Float64, bt2::Float64, pf1::Float64, pf2::Float64, btw::Float64, pfw::Float64}}(false, 500, (nrf = 0.4, bt1 = 0.2, bt2 = 0.6, pf1 = 0.4, pf2 = 0.5, btw = 0.6, pfw = 1.2), 5.0e-10, Bool[0 0 0], [0 0 0], [-1.0]), true)
If using a compound algorithm with a multi-point solve, then each of the solvers needs to have the multiple converged
and iterations
fields for each operating point, and the overall solve type needs to have a converged
and iterations
field for each solver and each operating point.
options = DuctAPE.set_options(;
solver_options=DuctAPE.ChainSolverOptions(;
solvers=[ # vector of solvers to use in poly-algorithm
DuctAPE.NLsolveOptions(;
algorithm=:anderson,
atol=1e-12,
iteration_limit=200,
converged=fill(false, (1, nop)), # flags for each operating point
iterations=zeros(Int, (1, nop)), # counters for each operating point
),
DuctAPE.MinpackOptions(;
atol=1e-12,
iteration_limit=100,
converged=fill(false, (1, nop)),
iterations=zeros(Int, (1, nop)),
),
],
converged=fill(false, (2, nop)), # flags for each solver and each operating point
iterations=zeros(Int, (2, nop)), # counts for each solver and each operating point
),
)
Options{Bool, DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}, Vector{Bool}, Float64, Vector{Int64}, Vector{String}, Vector{String}, DuctAPE.var"#56#61", IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}, ChainSolverOptions{Bool, Int64, Float64, DuctAPE.ExternalSolverOptions}, GridSolverOptions{Bool, Float64, Int64, Symbol}}(false, true, [1], DuctAPE.var"#56#61"(), true, false, 0.05, 1.0e-15, 2.220446049250313e-15, IntegrationOptions{GaussLegendre{Vector{Float64}, Vector{Float64}}, GaussLegendre{Vector{Float64}, Vector{Float64}}}(GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838]), GaussLegendre{Vector{Float64}, Vector{Float64}}([0.019855071751231856, 0.10166676129318664, 0.2372337950418355, 0.4082826787521751, 0.591717321247825, 0.7627662049581645, 0.8983332387068134, 0.9801449282487682], [0.05061426814518838, 0.11119051722668723, 0.15685332293894372, 0.18134189168918097, 0.18134189168918097, 0.15685332293894372, 0.11119051722668723, 0.05061426814518838])), DuctAPE.HeadsBoundaryLayerOptions{Bool, Float64, Float64, typeof(RK2), Int64, Float64, Float64, Float64, DuctAPE.RK, Nothing, Nothing}(false, true, true, false, 5000, 1.0e-6, nothing, nothing, 0.001, DuctAPE.RK(), DuctAPE.RK2, 3.0, 10, 10, 0.2, 0.2, false, 0.0, 0.0001, 0.0, false), Bool[0], ["outputs.jl"], false, ["outs"], GridSolverOptions{Bool, Float64, Int64, Symbol}(30, 3.0e-10, :newton, :forward, false, 3, Bool[0], [0], [0.0]), ChainSolverOptions{Bool, Int64, Float64, DuctAPE.ExternalSolverOptions}(DuctAPE.ExternalSolverOptions[NLsolveOptions{Bool, Float64, Int64, UnionAll, @NamedTuple{}, Symbol}(:anderson, 1.0e-12, 200, LineSearches.MoreThuente, NamedTuple(), Bool[0 0 0], [0 0 0], [-1.0]), MinpackOptions{Bool, Float64, Int64, Symbol}(:hybr, 1.0e-12, 100, Bool[0 0 0], [0 0 0], [-1.0])], Bool[0 0 0; 0 0 0], [0 0 0; 0 0 0], [-1.0]), true)