Analysis
DuctAPE.analyze_multipoint
— Functionanalyze_multipoint(
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;
solve_container_caching=nothing,
return_inputs=false,
)
Identical to the single analyze function assuming setup_analysis
has been called; except here we are running a single operating point for a multipoint analysis, and overwriting the operating point in the ducted_rotor with the explicit operating point input.
Process
DuctAPE.process
— Functionprocess(
solver_options::SolverOptionsType,
solve_parameter_cache_vector,
solve_parameter_cache_dims,
airfoils,
A_bb_LU,
solve_container_caching,
idmaps,
options,
)
Process (the step between pre-process and post-process) the solution, in other words: call the solver(s).
Arguments
solver_options::SolverOptionsType
: the solver options contained in the options object, used for dispatch.solve_parameter_cache_vector::Vector{Float}
: The vector cache for parameters used in the solve.solve_parameter_cache_dims::NamedTuple
: A named tuple containing the dimensions of the solve parameters.airfoils::NamedTuple
: The airfoils to be interpolated that are associated with each blade elementA_bb_LU::LinearAlgebra.LU
: The LU decomposition of the panel method LHS matrixsolve_container_caching::NamedTuple
: A named tuple containing the cache and dimensions for the intermediate solve values.idmaps::NamedTuple
: The set of index maps used in various solve sub-functionsoptions::Options
: User options
Returns
converged_states::Vector{Float}
: The output of a call toImplicitAD.implicit
DuctAPE.solve
— Functionsolve(sensitivity_parameters, const_cache; initial_guess=nothing)
A compact dispatch of solve
that automatically dispatches based on the solveroptions contained in constcache.
solve(
solver_options::SolverOptionsType,
sensitivity_parameters,
const_cache;
initial_guess=nothing,
)
Converge the residual, solving for the state variables that do so.
Arguments
solver_options::SolverOptionsType
: SolverOptionsType used for dispatchsensitivity_parameters::Vector{Float}
: Sensitivity parameters for solve (parameters passed in through ImplicitAD)const_cache::NamedTuple
: A named tuple containing constants and caching helpers.
Keyword Arguments
initial_guess=nothing::Vector{Float}
: An optional manually provided initial guess (contained in the sensitivity parameters anyway).
Returns
converged_states::Vector{Float}
: the states for which the residual has converged.
Internal Solvers
ModCSOR
DuctAPE.mod_COR_solver
— Functionmod_COR_solver(
r_fun!,
states,
B,
state_dims;
convergence_tolerance=1e-10,
iteration_limit=500,
relaxation_parameters=(; nrf=0.4, bt1=0.2, bt2=0.6, pf1=0.4, pf2=0.5, btw=0.6, pfw=1.2),
)
Modified DFDC-like CSOR solver that updates all states before relaxing Gamr and gamw.
Arguments:
r_fun!::function handle
: the residual function for the solver to useinitial_states::Vector{Float}
: the inital guess for the statesB::Vector{Float}
: number of blades on each rotor (used in state relaxation)state_dims::NamedTuple
: dimensions of the states (used in state relaxation)
Keyword Arguments
convergence_tolerance::type=1e-10
: absolute convergence toleranceiteration_limit::type=500
: maximum number of iterationsrelaxation_parameters::type=(; nrf=0.4, bt1=0.2, bt2=0.6, pf1=0.4, pf2=0.5, btw=0.6, pfw=1.2)
: parameters used in state relaxation
DuctAPE.relax_Gamr_mod!
— Functionrelax_Gamr_mod!(
Gamr,
r_Gamr_current,
r_Gamr_previous,
B;
nrf=0.4,
bt1=0.2,
bt2=0.6,
pf1=0.4,
pf2=0.5,
test=false,
)
Apply relaxed step to Gamr.
Arguments
Gamr::Array{Float}
: Array of rotor circulations (columns = rotors, rows = blade elements), updated in placeB::Vector{Float}
: number of blades on each rotorr_Gamr_current::Array{Float}
: Array of current iteration's differences in circulation valuesr_Gamr_previous::Array{Float}
: Array of previous iteration's differences in circulation values, updated in place
Keyword Arguments:
nrf::Float=0.4
: nominal relaxation factorbt1::Float=0.2
: backtrack factor 1bt2::Float=0.6
: backtrack factor 2pf1::Float=0.4
: press forward factor 1pf2::Float=0.5
: press forward factor 2
DuctAPE.relax_gamw_mod!
— Functionrelax_gamw_mod!(gamw, r_gamw_current, r_gamw_previous; nrf=0.4, btw=0.6, pfw=1.2)
Apply relaxed step to gamw.
Arguments
gamw::Array{Float}
: Array of rotor circulations (columns = rotors, rows = blade elements), updated in placer_gamw_current::Array{Float}
: Array of current iteration's differences in circulation valuesr_gamw_previous::Array{Float}
: Array of previous iteration's differences in circulation values, updated in place
Keyword Arguments
nrf::Float=0.4
: nominal relaxation factorbt1::Float=0.2
: backtrack factor 1bt2::Float=0.6
: backtrack factor 2pf1::Float=0.4
: press forward factor 1pf2::Float=0.5
: press forward factor 2
DuctAPE.update_states!
— Functionupdate_states!(states, r_current, r_previous, B, relaxation_parameters, state_dims)
Update states using DFDC-like relaxation methods.
Arguments
states::Vector{Float}
: current iteration states to updater_current::Vector{Float}
: current iteration residual valuesr_previous::Vector{Float}
: previous iteration residual valuesB::Vector{Float}
: number of blades for each rotorrelaxation_parameters::NamedTuple
: relaxation parametersstate_dims::NamedTuple
: dimensions of the state variablessolver_options::SolverOptionsType
: used for dispatch
CSOR
DuctAPE.apply_relaxation_schedule
— Functionapply_relaxation_schedule(
resid::AbstractVector, solver_options::TS
) where {TS<:SolverOptionsType}
Apply custom relaxation schedule to all relaxation factor inputs based on residual values.
Arguments
resid::AbstractVector{Float}
: current residual valuessolver_options::SolverOptionsType
: SolverOptions containing relaxation schedule
Returns
nrf::Float
: nominal relaxation factorbt1::Float
: backtrack factor 1bt2::Float
: backtrack factor 2pf1::Float
: press forward factor 1pf2::Float
: press forward factor 2
apply_relaxation_schedule(resid, nominal, schedule)
Apply custom relaxation schedule to a single relaxation factor input.
Arguments
resid::Float
: residual valuenominal::Float
: nominal relaxation valueschedule::AbstractVector{AbstractVector{Float}}
: values between which to interpolate to scale the nominal relaxation value.
Returns
rf::Float
: the updated relaxation factor
DuctAPE.update_CSOR_residual_values!
— Functionupdate_CSOR_residual_values!(
convergence_type::ConvergenceType, resid, maxBGamr, maxdeltaBGamr, maxdeltagamw, Vconv
)
Update CSOR residual values in place.
Arguments
convergence_type::ConvergenceType
: used for dispatch of relative or absolute residual values.resid::Vector{Float}
: residual values modified in placemaxBGamr::Float
: Maximum value of B*Gamr among all blade elementsmaxdeltaBGamr::Float
: Maximum change in B*Gamr between iterations among all blade elementsmaxdeltagamw::Vector{Float}
: Maximum change in gamw among all wake nodes (one element)Vconv::Float
: Reference velocity upon which the relative convergence criteria is based (one element)
DuctAPE.check_CSOR_convergence!
— Functioncheck_CSOR_convergence!(
conv, resid; f_circ=1e-3, f_dgamw=2e-4, convergence_type=Relative(), verbose=false
)
Description
Arguments
conv::Vector{Float}
: container holding convergence flagresid::Vector{Float}
: residual vector
Keyword Arguments
f_circ::Float=1e-3
: convergence criteria for circulation residualf_dgamw::Float=2e-4
: convergence criteria for wake strength residualconvergence_type::ConvergenceType=Relative()
: convergence type (absolute or relative) for print statementsverbose::Bool=false
: flag for verbose print statements
DuctAPE.relax_Gamr!
— Functionrelax_Gamr!(
Gamr,
delta_prev_mat,
delta_mat,
maxBGamr,
maxdeltaBGamr,
B;
nrf=0.4,
bt1=0.2,
bt2=0.6,
pf1=0.4,
pf2=0.5,
test=false,
)
Apply relaxed step to Gamr.
Arguments
Gamr::Array{Float}
: Array of rotor circulations (columns = rotors, rows = blade elements), updated in placedelta_prev_mat::Array{Float}
: Array of previous iteration's differences in circulation values, updated in placedelta_mat::Array{Float}
: Array of current iteration's differences in circulation valuesmaxBGamr::Array{Float}
: stores value of maximum B*Gamr for each rotormaxdeltaBGamr::Array{Float}
: stores value of maximum change in B*Gamr for each rotorB::Vector{Float}
: number of blades on each rotornrf::Float=0.4
: nominal relaxation factorbt1::Float=0.2
: backtrack factor 1bt2::Float=0.6
: backtrack factor 2pf1::Float=0.4
: press forward factor 1pf2::Float=0.5
: press forward factor 2
DuctAPE.relax_gamw!
— Functionrelax_gamw!(
gamw, delta_prev, delta, maxdeltagamw; nrf=0.4, btw=0.6, pfw=1.2, test=false
)
Apply relaxed step to gamw.
Arguments
gamw::Array{Float}
: Array of rotor circulations (columns = rotors, rows = blade elements), updated in placedelta_prev_mat::Array{Float}
: Array of previous iteration's differences in circulation values, updated in placedelta_mat::Array{Float}
: Array of current iteration's differences in circulation valuesmaxdeltagamw::Array{Float}
: Single element array that gets updated with the new maximum change in gamw.
Keyword Arguments:
nrf::Float=0.4
: nominal relaxation factorbt1::Float=0.2
: backtrack factor 1bt2::Float=0.6
: backtrack factor 2pf1::Float=0.4
: press forward factor 1pf2::Float=0.5
: press forward factor 2
Residuals
ModCSOR
DuctAPE.mod_CSOR_residual!
— Functionmod_CSOR_residual!(r, current_states, inputs, constants)
Modified DFDC-like CSOR residual that does not include any relaxation within the residual calculation.
Arguments
r::Vector{Float}
: solve residualcurrent_states::Vector{Float}
: solve statesinputs::Vector{Float}
: solve inputsconstants::Vector{Float}
: solve constants
DuctAPE.estimate_CSOR_states!
— Functionestimate_CSOR_states!(
solve_containers,
Gamr,
sigr,
gamw,
operating_point,
ivr,
ivw,
linsys,
blade_elements,
wakeK,
idmaps;
verbose=false,
)
Estimate states for modified CSOR solver.
Arguments:
solve_containers::NamedTuple
: cache for intermediate solve valuesGamr::type
: Blade element circulation strengthssigr::type
: Rotor source panel strengthsgamw::type
: Wake vortex panel strengthsoperating_point::NamedTuple
: Named tuple containing operating_point informationivr::NamedTuple
: unit induced velocities on rotor(s)ivw::NamedTuple
: unit induced velocities on wakelinsys::NamedTuple
: vectors and matricies comprising the panel method linear systemblade_elements::NamedTuple
: blade element geometry and airfoil polar informationwakeK::Vector{Float}
: geometric constants used in caculating wake strengthsidmaps::NamedTuple
: index maps used throughout solve
CSOR
DuctAPE.CSOR_residual!
— FunctionCSOR_residual!(resid, state_variables, sensitivity_parameters, constants)
The in-place residual used for the CSOR solve method.
Arguments
resid::Vector{Float}
: In-place residual.state_variables::Vector{Float}
: The state variablessensitivity_parameters::Vector{Float}
: The parameters to which the solution is sensitive.constants::NamedTuple
: Various constants required in the solve
Returns
state_variables::Vector{Float}
: The state variables (modified in place)
DuctAPE.compute_CSOR_residual!
— Functioncompute_CSOR_residual!(
resid,
solver_options,
solve_containers,
Gamr,
sigr,
gamw,
operating_point,
ivr,
ivw,
linsys,
blade_elements,
wakeK,
idmaps;
verbose=false,
)
Description
Arguments
resid::Vector{Float}
: the residual vectorsolver_options::SolverOptionsType
: solver options (used for convergence criteria)solve_containers::NamedTuple
: cache for intermediate solve valuesGamr::type
: Blade element circulation strengthssigr::type
: Rotor source panel strengthsgamw::type
: Wake vortex panel strengthsoperating_point::NamedTuple
: Named tuple containing operating_point informationivr::NamedTuple
: unit induced velocities on rotor(s)ivw::NamedTuple
: unit induced velocities on wakelinsys::NamedTuple
: vectors and matricies comprising the panel method linear systemblade_elements::NamedTuple
: blade element geometry and airfoil polar informationwakeK::Vector{Float}
: geometric constants used in caculating wake strengthsidmaps::NamedTuple
: index maps used throughout solve
Keyword Arguments
verbose::Bool=false
: Flag to print verbose statements
External Solvers
DuctAPE.system_residual
— Functionsystem_residual(state_variables, sensitivity_parameters, constants)
The residual function for external solvers.
Arguments
state_variables::Vector{Float}
: the state variablessensitivity_parameters::Vector{Float}
: parameters to which the solution derivatives are sensitiveconstants::NamedTuple
: parameters to which the solution derivatives are constant
Returs
resid::Vector{Float}
: residual vector
DuctAPE.system_residual!
— Functionsystem_residual!(resid, state_variables, sensitivity_parameters, constants)
In-place version of system_residual.
DuctAPE.update_system_residual!
— Functionupdate_system_residual!(
solver_options::SolverOptionsType
resid,
vz_est,
vz_rotor,
vtheta_est,
vtheta_rotor,
Cm_est,
Cm_wake,
solve_parameter_cache_dims,
)
Update the residual for external solvers.
Arguments
- `solver_options::SolverOptionsType
resid::Vector{Float}
: residual vectorvz_est::Vector{Float}
: axial induced rotor velocity estimate containervz_rotor::Vector{Float}
: axial induced rotor velocity state containervtheta_est::Vector{Float}
: tangential induced rotor velocity estimate containervtheta_rotor::Vector{Float}
: tangential induced rotor velocity state containerCm_est::Vector{Float}
: absolute meridional wake control point velocity estimate containerCm_wake::Vector{Float}
: absolute meridional wake control point velocity state containersolve_parameter_cache_dims::Vector{Float}
: dimensions of state vectors to use in accessing the residual vector
DuctAPE.estimate_states!
— Functionestimate_states!(
solve_containers,
vz_rotor,
vtheta_rotor,
Cm_wake,
operating_point,
ivr,
ivw,
linsys,
blade_elements,
wakeK,
idmaps;
verbose=false,
)
Estimate velocity states.
Arguments
solve_containers::NamedTuple
: cache for intermediate values in solvevz_rotor::Vector{Float}
: axial induced rotor velocity state containervtheta_rotor::Vector{Float}
: tangential induced rotor velocity state containerCm_wake::Vector{Float}
: absolute meridional wake control point velocity state containeroperating_point::NamedTuple
: Named tuple containing operating_point informationivr::NamedTuple
: unit induced velocities on rotor(s)ivw::NamedTuple
: unit induced velocities on wakelinsys::NamedTuple
: vectors and matricies comprising the panel method linear systemblade_elements::NamedTuple
: blade element geometry and airfoil polar informationwakeK::Vector{Float}
: geometric constants used in caculating wake strengthsidmaps::NamedTuple
: index maps used throughout solve
Keyword Arguments
verbose::Bool=false
: flag for verbose print statements
Solve Utilities
DuctAPE.extract_initial_guess
— Functionextract_initial_guess(
solver_options::SolverOptionsType, sensitivity_parameters, state_dims
)
Extract initial guess from the solve parameters cache vector.
Arguments
solver_options::SolverOptionsType
: used for dispatchsensitivity_parameters::Vector{Float}
: vector form of solve parameter cache passed into the solver.state_dims::NamedTuple
: dimensions and indices of state variables within the solve parameter cache vector
Returns
- initial_guess::Vector{Float}` : a vector of the solver initial guess
DuctAPE.extract_state_variables
— Functionextract_state_variables(solver_options::SolverOptionsType, vars, dims)
Reshape the state variables from a single vector, to multiple arrays.
Arguments
Returns if solver_options <: CSORSolverOptions
Gamr::type
: Blade element circulation strengthssigr::type
: Rotor source panel strengthsgamw::type
: Wake vortex panel strengths
Returns if solver_options <: Union{ExternalSolverOptions, PolyAlgorithmOptions}
vz_rotor::Vector{Float}
: axial induced rotor velocity state containervtheta_rotor::Vector{Float}
: tangential induced rotor velocity state containerCm_wake::Vector{Float}
: absolute meridional wake control point velocity state container