Module: lambertSolver
Executive Summary
This module solves Lambert’s problem using either the algorithm by R.H. Gooding
(longer, extensive report available here ) or the
algorithm by D. Izzo. Given two position vectors
The algorithms by Gooding and Izzo provide solutions for elliptic, parabolic and hyperbolic transfer orbits. In the zero-revolution case, exactly one solution exists. In the multi-revolution case (meaning that the transfer orbit completes at least one complete revolution about the central body), either two or zero solutions exist, depending on whether the requested transfer time is greater or less than the minimum time-of-flight for the given number of revolutions. The LambertSolutionMsgPayload output message includes two solutions. If a solution does not exist, the corresponding velocity vectors in the message are equal to the zero vector, and the “validity” flag is equal to zero. Otherwise, the “validity” flag is set to 1.
Message Connection Descriptions
The following table lists all the module input and output messages. The module msg connection is set by the user from python. The msg type contains a link to the message structure definition, while the description provides information on what this message is used for.
Msg Variable Name |
Msg Type |
Description |
---|---|---|
lambertProblemInMsg |
lambert problem setup input message |
|
lambertSolutionOutMsg |
lambert problem solution output message |
|
lambertPerformanceOutMsg |
lambert problem performance message (additional information about the solution process) |
Module Assumptions and Limitations
The algorithms only compute solutions for a positive time-of-flight, and for positive transfer angles (meaning that the
true anomaly of
An edge case exists for a transfer angle of 0 or 180 degrees, as the two position vectors do not define a plane, so an infinite number of solutions exist. The module checks if the angle between the two position vectors is smaller than the threshold “alignmentThreshold”. In this case, the computed velocity vectors are set equal to the zero vector and the validity flag of the solution is set to zero.
While the module also works for parabolic transfer orbits, the solutions for orbits that are very close to - but not exactly - parabolic, may not be as accurate.
User Guide
The module is first initialized as follows:
lambertModule = lambertSolver.LambertSolver()
lambertModule.ModelTag = "lambertSolver"
lambertModule.setAlignmentThreshold(1.0) # module defaults this value to 1.0 degrees if not specified
unitTestSim.AddModelToTask(unitTaskName, lambertModule)
The lambert problem input message is either created as a standalone message in python
lambertProblemInMsgData = messaging.LambertProblemMsgPayload()
lambertProblemInMsgData.solverMethod = messaging.IZZO
lambertProblemInMsgData.r1_N = np.array([10000. * 1000, 0. ,0.])
lambertProblemInMsgData.r2_N = np.array([0., 8000. * 1000,0.])
lambertProblemInMsgData.transferTime = 10000.
lambertProblemInMsgData.mu = 3.986004418e14
lambertProblemInMsgData.numRevolutions = 0
lambertProblemInMsg = messaging.LambertProblemMsg().write(lambertProblemInMsgData)
or obtained from another FSW module. The lambert problem input message is then connected.
lambertModule.lambertProblemInMsg.subscribeTo(lambertProblemInMsg)
-
class LambertSolver : public SysModel
- #include <lambertSolver.h>
This module solves Lambert’s problem using either the Gooding or the Izzo algorithm.
Public Functions
-
LambertSolver()
This is the constructor for the module class. It sets default variable values and initializes the various parts of the model
-
~LambertSolver()
Module Destructor
-
void Reset(uint64_t currentSimNanos) override
This method is used to reset the module and checks that required input messages are connected.
- Parameters:
currentSimNanos – current simulation time in nano-seconds
- Returns:
void
-
void UpdateState(uint64_t currentSimNanos) override
This is the main method that gets called every time the module is updated. It computes the solution of Lambert’s problem.
- Parameters:
currentSimNanos – current simulation time in nano-seconds
- Returns:
void
-
void setAlignmentThreshold(const double value)
setter for
alignmentThreshold
-
inline double getlignmentThreshold() const
getter for
alignmentThreshold
Public Members
-
ReadFunctor<LambertProblemMsgPayload> lambertProblemInMsg
lambert problem input message
-
Message<LambertSolutionMsgPayload> lambertSolutionOutMsg
lambert solution output message
-
Message<LambertPerformanceMsgPayload> lambertPerformanceOutMsg
lambert performance output message
-
BSKLogger bskLogger
BSK Logging.
Private Functions
-
void readMessages()
This method reads the input messages each call of updateState. It also checks if the message contents are valid for this module.
- Returns:
void
-
void writeMessages(uint64_t currentSimNanos)
This method writes the output messages each call of updateState
- Parameters:
currentSimNanos – current simulation time in nano-seconds
- Returns:
void
-
void problemGeometry()
This method computes the problem geometry for the given parameters of Lambert’s problem. The orbit frame is also determined.
- Returns:
void
-
std::array<double, 2> goodingInitialGuess(double lambda, double T)
This method computes the initial guess for the free variable x using the Gooding procedure
- Parameters:
lam – lambda parameter that defines the problem geometry
T – non-dimensional time-of-flight
- Returns:
std::array<double, 2>
-
std::array<double, 2> izzoInitialGuess(double lambda, double T)
This method computes the initial guess for the free variable x using the Izzo procedure
- Parameters:
lam – lambda parameter that defines the problem geometry
T – non-dimensional time-of-flight
- Returns:
std::array<double, 2>
-
void findx()
This method finds the free variable x that satisfies the requested time of flight TOF.
- Returns:
void
-
std::array<Eigen::Vector3d, 2> computeVelocities(double x)
This method computes the velocities at the initial and final position for a given free variable x
- Parameters:
x – free variable of Lambert’s problem that satisfies the given time of flight
- Returns:
std::array<Eigen::Vector3d, 2>
-
double x2tof(double x, int N, double lam)
This method computes the non-dimensional time of flight (TOF) for a given x
- Parameters:
x – free variable of Lambert’s problem that satisfies the given time of flight
N – number of revolutions
lam – lambda parameter that defines the problem geometry
- Returns:
double
-
std::array<double, 3> dTdx(double x, double T, double lam)
This method computes the derivatives of the time of flight curve T(x)
- Parameters:
x – free variable of Lambert’s problem that satisfies the given time of flight
T – requested non-dimensional time-of-flight
lam – lambda parameter that defines the problem geometry
- Returns:
std::array<double, 3>
-
std::array<double, 3> householder(double T, double x0, int N)
This method includes a 3rd order householder root-finder to find the free variable x that satisfies the time-of-flight constraint.
- Parameters:
T – requested non-dimensional time-of-flight
x0 – initial guess for x free variable of Lambert’s problem that satisfies the given time of flight
N – number of revolutions
- Returns:
std::array<double, 3>
-
std::array<double, 3> halley(double T, double x0, int N)
This method includes a halley root-finder (2nd order householder) to find the free variable x that satisfies the time-of-flight constraint
- Parameters:
T – requested non-dimensional time-of-flight
x0 – initial guess for x free variable of Lambert’s problem that satisfies the given time of flight
N – number of revolutions
- Returns:
std::array<double, 3>
-
double getTmin(double T0M, int N)
This method computes the minimum non-dimensional time-of-flight Tmin such that solutions exist for the multi-revolution case using a halley root-finder
- Parameters:
T0M – initial guess for Tmin
N – number of revolutions
- Returns:
double
-
double hypergeometricF(double z)
This method computes the hypergeometric function 2F1(a,b,c,z)
- Parameters:
z – argument of hypergeometric function
- Returns:
double
Private Members
-
double alignmentThreshold = 1.0
[deg] minimum angle between position vectors so they are not too aligned.
-
SolverMethod solverMethod
lambert solver algorithm (GOODING or IZZO)
-
Eigen::Vector3d r1_N
position vector at t0
-
Eigen::Vector3d r2_N
position vector at t1
-
double transferTime = {}
time of flight between r1_N and r2_N (t1-t0)
-
double mu = {}
gravitational parameter
-
int numberOfRevolutions = {}
number of revolutions
-
double TOF = {}
non-dimensional time-of-flight constraint
-
double lambda = {}
parameter of Lambert”s problem that defines problem geometry
-
bool multiRevSolution = {}
boolean flag if multi-revolution solutions exist or not
-
bool noSolution = {}
boolean flag if no solution should be returned (in case of 180 deg transfer angle)
-
std::array<Eigen::Vector3d, 3> Oframe1
array containing the orbit frame unit vectors at t0
-
std::array<Eigen::Vector3d, 3> Oframe2
array containing the orbit frame unit vectors at t1
-
std::array<Eigen::Vector3d, 2> vvecs
array containing the velocity vector solutions at t0 and t1
-
std::array<Eigen::Vector3d, 2> vvecsSol2
array containing the velocity vector solutions at t0 and t1 (sol 2)
-
double X = {}
solution for free variable of Lambert’s problem
-
double XSol2 = {}
second solution for free variable of Lambert’s problem
-
int numIter = {}
number of root finder iterations to find X
-
int numIterSol2 = {}
number of root finder iterations to find X_sol2
-
double errX = {}
difference in X between last and second-to-last iteration
-
double errXSol2 = {}
difference in X between last and second-to-last iteration (for X_sol2)
-
LambertSolver()