return GSL_SUCCESS; } struct Parameters{ }; TEST_CASE( "Print Cartesian state to Two-Line-Elements table header", "[print]" ) { // Set expected output string. std::ostringstream tableHeader; tableHeader << "# a e i AoP RAAN " << "TA f1 f2 f3 f4 " << "f5 f6 " << std::endl; REQUIRE( printCartesianToTleSolverStateTableHeader( ) == tableHeader.str( ) ); } TEST_CASE( "Print Cartesian state to Two-Line-Elements solver state", "[print]" ) { // Set initial guess for GSL solver. gsl_vector* initialGuess = gsl_vector_alloc( 6 ); gsl_vector_set( initialGuess, 0, 0.1 ); gsl_vector_set( initialGuess, 1, 0.2 ); gsl_vector_set( initialGuess, 2, 0.3 ); gsl_vector_set( initialGuess, 3, 0.4 ); gsl_vector_set( initialGuess, 4, 0.5 ); gsl_vector_set( initialGuess, 5, 0.6 ); // Set up dummy GSL solver. Parameters parameters;
const Tle convertCartesianStateToTwoLineElements( const Vector6& cartesianState, const DateTime& epoch, std::string& solverStatusSummary, int& numberOfIterations, const Tle& referenceTle, const Real earthGravitationalParameter, const Real earthMeanRadius, const Real absoluteTolerance, const Real relativeTolerance, const int maximumIterations ) { // Store reference TLE as the template TLE and update epoch. Tle templateTle = referenceTle; templateTle.updateEpoch( epoch ); // Set up parameters for residual function. CartesianToTwoLineElementsParameters< Vector6 > parameters( cartesianState, templateTle ); // Set up residual function. gsl_multiroot_function cartesianToTwoLineElementsFunction = { &computeCartesianToTwoLineElementResiduals< Real, Vector6 >, 6, ¶meters }; // Compute current state in Keplerian elements, for use as initial guess for the TLE mean // elements. const Vector6 initialKeplerianElements = astro::convertCartesianToKeplerianElements( parameters.targetState, earthGravitationalParameter ); // Compute initial guess for TLE mean elements. const Vector6 initialTleMeanElements = computeInitialGuessTleMeanElements( initialKeplerianElements, earthGravitationalParameter ); // Set initial guess. gsl_vector* initialGuessTleMeanElements = gsl_vector_alloc( 6 ); for ( int i = 0; i < 6; i++ ) { gsl_vector_set( initialGuessTleMeanElements, i, initialTleMeanElements[ i ] ); } // Set up solver type (derivative free). const gsl_multiroot_fsolver_type* solverType = gsl_multiroot_fsolver_hybrids; // Allocate memory for solver. gsl_multiroot_fsolver* solver = gsl_multiroot_fsolver_alloc( solverType, 6 ); // Set solver to use residual function with initial guess for TLE mean elements. gsl_multiroot_fsolver_set( solver, &cartesianToTwoLineElementsFunction, initialGuessTleMeanElements ); // Declare current solver status and iteration counter. int solverStatus = false; int counter = 0; // Set up buffer to store solver status summary table. std::ostringstream summary; // Print header for summary table to buffer. summary << printCartesianToTleSolverStateTableHeader( ); do { // Print current state of solver for summary table. summary << printCartesianToTleSolverState( counter, solver ); // Increment iteration counter. ++counter; // Execute solver iteration. solverStatus = gsl_multiroot_fsolver_iterate( solver ); // Check if solver is stuck; if it is stuck, break from loop. if ( solverStatus ) { solverStatusSummary = summary.str( ); throw std::runtime_error( "ERROR: Non-linear solver is stuck!" ); } // Check if root has been found (within tolerance). solverStatus = gsl_multiroot_test_delta( solver->dx, solver->x, absoluteTolerance, relativeTolerance ); } while ( solverStatus == GSL_CONTINUE && counter < maximumIterations ); // Save number of iterations. numberOfIterations = counter - 1; // Print final status of solver to buffer. summary << std::endl; summary << "Status of non-linear solver: " << gsl_strerror( solverStatus ) << std::endl; summary << std::endl; // Write buffer contents to solver status summary string. solverStatusSummary = summary.str( ); // Generate TLE with converged mean elements. Tle virtualTle = templateTle; Real convergedMeanEccentricity = gsl_vector_get( solver->x, 2 ); if ( convergedMeanEccentricity < 0.0 ) { convergedMeanEccentricity = std::fabs( gsl_vector_get( solver->x, 2 ) ); } if ( convergedMeanEccentricity > 0.999 ) { convergedMeanEccentricity = 0.99; } virtualTle.updateMeanElements( sml::computeModulo( std::fabs( gsl_vector_get( solver->x, 0 ) ), 180.0 ), sml::computeModulo( gsl_vector_get( solver->x, 1 ), 360.0 ), convergedMeanEccentricity, sml::computeModulo( gsl_vector_get( solver->x, 3 ), 360.0 ), sml::computeModulo( gsl_vector_get( solver->x, 4 ), 360.0 ), std::fabs( gsl_vector_get( solver->x, 5 ) ) ); // Free up memory. gsl_multiroot_fsolver_free( solver ); gsl_vector_free( initialGuessTleMeanElements ); return virtualTle; }