Periodic Orbit Computation
CRTBPNaturalMotion.GeneralPeriodicOrbit — Type
GeneralPeriodicOrbitA general representation of a periodic orbit.
Fields
x0::SVector{6, Float64}: The initial conditions for the periodic orbit.P::Float64: The period of the periodic orbit.S::Float64: The arc length along the periodic orbit.mu::Float64: The gravitational parameter of the CRTBP.
CRTBPNaturalMotion.TypeAPeriodicOrbit — Type
TypeAPeriodicOrbit(
rx, rz, vy, mu[, TU, DU];
N_cross = 1,
constraint = (:x_start_coordinate, 0.8),
ode_solver = Vern9(),
ode_reltol = 1e-12,
ode_abstol = 1e-12,
nl_solver = SimpleTrustRegion(autodiff = nothing, nlsolve_update_rule = Val(true)),
)CRTBPNaturalMotion.TypeAPeriodicOrbit — Type
TypeAPeriodicOrbitA representation of a CRTBP periodic orbit of Type A.
Fields
u0::SVector{3, Float64}: The initial conditions for the periodic orbit.P::Float64: The period of the periodic orbit.S::Float64: The arc length along the periodic orbit.N_cross::Int: The number of crossings through x-z plane in half period.mu::Float64: The gravitational parameter of the CRTBP.
CRTBPNaturalMotion.TypeAPeriodicOrbit — Method
TypeAPeriodicOrbit(
orbit::TypeAPeriodicOrbit;
constraint = (:x_start_coordinate, 0.8),
ode_solver = Vern9(),
ode_reltol = 1e-12,
ode_abstol = 1e-12,
nl_solver = SimpleTrustRegion(autodiff = nothing, nlsolve_update_rule = Val(true)),
)CRTBPNaturalMotion.generate_periodic_orbit_cheb_approximation — Method
generate_periodic_orbit_cheb_approximation(
orbit::AbstractPeriodicOrbit, τ1_npoints::Int, τ1_order::Int, PT::ParameterizationType;
solver = Vern9(),
reltol = 1e-14,
abstol = 1e-14,
) -> FastChebInterpolation{FastChebInterp.ChebPoly{1, SVector{6, Float64}, Float64}}Generate a Chebyshev approximation of the periodic orbit in terms of PT employing τ1_npoints fit to a Chebyshev polynomial of order τ1_order.
Arguments
orbit::AbstractPeriodicOrbit: The periodic orbit.τ1_npoints::Int: The number of points to fit the Chebyshev polynomial.τ1_order::Int: The order of the Chebyshev polynomial.PT::ParameterizationType: The type of parameterization to use (TimeorArcLength).
Keyword Arguments
solver::OrdinaryDiffEq.AbstractODESolver: The ODE solver to use for the integration.reltol::Real: The relative tolerance for the ODE solver.abstol::Real: The absolute tolerance for the ODE solver.
Returns
FastChebInterpolation{FastChebInterp.ChebPoly{1, SVector{6, Float64}, Float64}}: The Chebyshev approximation.
CRTBPNaturalMotion.generate_periodic_orbit_cheb_interpolation — Method
generate_periodic_orbit_cheb_interpolation(
orbit::AbstractPeriodicOrbit, τ1_order::Int, PT::ParameterizationType;
solver = Vern9(),
reltol = 1e-14,
abstol = 1e-14,
) -> FastChebInterpolation{FastChebInterp.ChebPoly{1, SVector{6, Float64}, Float64}}Generate a Chebyshev interpolation of the periodic orbit in terms of PT employing a Chebyshev polynomial of order τ1_order.
Arguments
orbit::AbstractPeriodicOrbit: The periodic orbit.τ1_order::Int: The order of the Chebyshev polynomial.PT::ParameterizationType: The type of parameterization to use (TimeorArcLength).
Keyword Arguments
solver::OrdinaryDiffEq.AbstractODESolver: The ODE solver to use for the integration.reltol::Real: The relative tolerance for the ODE solver.abstol::Real: The absolute tolerance for the ODE solver.
Returns
FastChebInterpolation{FastChebInterp.ChebPoly{1, SVector{6, Float64}, Float64}}: The Chebyshev interpolation.
CRTBPNaturalMotion.get_full_orbit — Method
get_full_orbit(orbit::AbstractPeriodicOrbit)Get the full orbit of the periodic orbit.
Arguments
orbit::AbstractPeriodicOrbit: The periodic orbit.
CRTBPNaturalMotion.correct_typeA_initial_conditions — Method
correct_typeA_initial_conditions(
rx_guess, rz_guess, vy_guess, mu;
N_cross = 1,
constraint = (:x_start_coordinate, 0.0),
ode_solver = Vern9(),
ode_reltol = 1e-14,
ode_abstol = 1e-14,
nl_solver = SimpleTrustRegion(autodiff = nothing, nlsolve_update_rule = Val(true)),
)Correct the initial conditions of a type-A periodic orbit using a forward shooting method.
Arguments
rx_guess::Real: Initial guess for the x-coordinate of the initial state.rz_guess::Real: Initial guess for the z-coordinate of the initial state.vy_guess::Real: Initial guess for the y-velocity of the initial state.mu::Real: Gravitational parameter of the CRTBP.
Keyword Arguments
N_cross::Int = 1: Number of crossings of the xz-plane to consider in the half period.constraint::Tuple{Symbol,Real}: Tuple of the constraint type and value corresponding to the final constraint considered to provide a fully defined system of nonlinear eqautions. The constraint type can be either:x_start_coordinateor:jacobi_integral. The constraint value is the value that the constraint should be equal to for the desired periodic orbit.ode_solver::OrdinaryDiffEq.AbstractODESolver: The ODE solver to use for the integration.ode_reltol::Real: The relative tolerance for the ODE solver.ode_abstol::Real: The absolute tolerance for the ODE solver.nl_solver::NonlinearSolve.AbstractNLsolveSolver: The nonlinear solver to use for the root-finding problem. Note, there seems to be a strange bug with the NonlinearSolve.jlTrustRegionandNewtonRaphsonsolvers when using out-of-place Jacobian defined to return an SMatrix. Therefore, it is recommened to use either the defaultSimpleTrustRegionsolver or theSimpleNewtonRaphsonsolver (which are specifically tailored towards StaticArrays.jl arrays).