Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Decouple KeplerianOrbit from GravityEffector #5

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
105 changes: 43 additions & 62 deletions src/architecture/utilities/keplerianOrbit.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -18,122 +18,107 @@
*/

#include "keplerianOrbit.h"
#include "architecture/utilities/astroConstants.h"
#include <architecture/utilities/avsEigenSupport.h>
#include <architecture/utilities/linearAlgebra.h>

/*! This constructor initialized to an arbitrary orbit */
KeplerianOrbit::KeplerianOrbit()
{
this->semi_major_axis = 1E5;
this->eccentricity = 1E-5;
this->inclination = 0.0;
this->true_anomaly = 0.0;
this->argument_of_periapsis = 0.0;
this->right_ascension = 0.0;
this->mu = MU_EARTH;
this->change_orbit();
return;
}

/*! The constructor requires orbital elements and a planet */
KeplerianOrbit::KeplerianOrbit(classicElements oe, GravBodyData* planet){
this->set_planet(planet);
this->semi_major_axis = oe.a;
this->eccentricity = oe.e;
this->inclination = oe.i;
this->true_anomaly = oe.f;
this->argument_of_periapsis = oe.omega;
this->right_ascension = oe.Omega;
this->mu = planet->mu;
/*! The constructor requires orbital elements and a gravitational constant value */
KeplerianOrbit::KeplerianOrbit(classicElements oe, const double mu) : mu(mu),
semi_major_axis(oe.a),
eccentricity(oe.e),
inclination(oe.i),
argument_of_periapsis(oe.omega),
right_ascension(oe.Omega),
true_anomaly(oe.f){
this->change_orbit();
return;
}

/*! The copy constructor works with python copy*/
KeplerianOrbit::KeplerianOrbit(const KeplerianOrbit &orig){
this->semi_major_axis = orig.semi_major_axis;
this->eccentricity = orig.eccentricity;
this->inclination = orig.inclination;
this->true_anomaly = orig.true_anomaly;
this->argument_of_periapsis = orig.argument_of_periapsis;
this->right_ascension = orig.right_ascension;
this->mu = orig.mu;
this->planet = orig.planet;
KeplerianOrbit::KeplerianOrbit(const KeplerianOrbit &orig) : mu(orig.mu),
semi_major_axis(orig.a()),
eccentricity(orig.e()),
inclination(orig.i()),
argument_of_periapsis(orig.omega()),
right_ascension(orig.RAAN()),
true_anomaly(orig.f())
{
this->change_orbit();
return;
}

/*! Generic Destructor */
KeplerianOrbit::~KeplerianOrbit()
{
return;
}

/*!
body position vector relative to planet
*/
Eigen::Vector3d KeplerianOrbit::r_BP_P(){
Eigen::Vector3d KeplerianOrbit::r_BP_P() const {
return this->position_BP_P;
}

/*!
body velocity vector relative to planet
*/
Eigen::Vector3d KeplerianOrbit::v_BP_P(){
Eigen::Vector3d KeplerianOrbit::v_BP_P() const{
return this->velocity_BP_P;

}

/*!
angular momentum of body relative to planet
*/
Eigen::Vector3d KeplerianOrbit::h_BP_P(){
Eigen::Vector3d KeplerianOrbit::h_BP_P() const{
return this->orbital_angular_momentum_P;
}

/*! return mean anomaly angle */
double KeplerianOrbit::M(){return this->mean_anomaly;}
double KeplerianOrbit::M() const {return this->mean_anomaly;}
/*! return mean orbit rate */
double KeplerianOrbit::n(){return this->mean_motion;}; //!< return mean orbit rate
double KeplerianOrbit::n() const {return this->mean_motion;}; //!< return mean orbit rate
/*! return orbit period */
double KeplerianOrbit::P(){return this->orbital_period;}; //!< return orbital period
double KeplerianOrbit::P() const {return this->orbital_period;}; //!< return orbital period
/*! return true anomaly */
double KeplerianOrbit::f(){return this->true_anomaly;}; //!< return true anomaly
double KeplerianOrbit::f() const {return this->true_anomaly;}; //!< return true anomaly
/*! return true anomaly rate */
double KeplerianOrbit::fDot(){return this->true_anomaly_rate;};
double KeplerianOrbit::fDot() const {return this->true_anomaly_rate;};
/*! return right ascencion of the ascending node */
double KeplerianOrbit::RAAN(){return this->right_ascension;};
double KeplerianOrbit::RAAN() const {return this->right_ascension;};
/*! return argument of periapses */
double KeplerianOrbit::omega(){return this->argument_of_periapsis;};
double KeplerianOrbit::omega() const {return this->argument_of_periapsis;};
/*! return inclination angle */
double KeplerianOrbit::i(){return this->inclination;};
double KeplerianOrbit::i() const {return this->inclination;};
/*! return eccentricty */
double KeplerianOrbit::e(){return this->eccentricity;};
double KeplerianOrbit::e() const {return this->eccentricity;};
/*! return semi-major axis */
double KeplerianOrbit::a(){return this->semi_major_axis;};
double KeplerianOrbit::a() const {return this->semi_major_axis;};
/*! return orbital angular momentum magnitude */
double KeplerianOrbit::h(){return this->h_BP_P().norm();};
double KeplerianOrbit::h() const {return this->h_BP_P().norm();};
/*! return orbital energy */
double KeplerianOrbit::Energy(){return this->orbital_energy;};
/*! return orbit radius */
double KeplerianOrbit::r(){return this->r_BP_P().norm();};
double KeplerianOrbit::r() const {return this->r_BP_P().norm();};
/*! return velocity magnitude */
double KeplerianOrbit::v(){return this->v_BP_P().norm();};
double KeplerianOrbit::v() const {return this->v_BP_P().norm();};
/*! return radius at apoapses */
double KeplerianOrbit::r_a(){return this->r_apogee;};
double KeplerianOrbit::r_a() const {return this->r_apogee;};
/*! return radius at periapses */
double KeplerianOrbit::r_p(){return this->r_perigee;};
double KeplerianOrbit::r_p() const {return this->r_perigee;};
/*! return flight path angle */
double KeplerianOrbit::fpa(){return this->flight_path_angle;};
double KeplerianOrbit::fpa() const {return this->flight_path_angle;};
/*! return eccentric anomaly angle */
double KeplerianOrbit::E(){return this->eccentric_anomaly;};
double KeplerianOrbit::E() const {return this->eccentric_anomaly;};
/*! return semi-latus rectum or the parameter */
double KeplerianOrbit::p(){return this->semi_parameter;};
double KeplerianOrbit::p() const {return this->semi_parameter;};
/*! return radius rate */
double KeplerianOrbit::rDot(){return this->radial_rate;};
double KeplerianOrbit::rDot() const {return this->radial_rate;};
/*! return escape velocity */
double KeplerianOrbit::c3(){return this->v_infinity;};
double KeplerianOrbit::c3() const {return this->v_infinity;};

/*! set semi-major axis */
void KeplerianOrbit::set_a(double a){this->semi_major_axis = a; this->change_orbit();};
Expand Down Expand Up @@ -179,7 +164,6 @@ void KeplerianOrbit::change_orbit(){
this->orbital_energy = -this->mu / 2 / this->a();
this->r_apogee = this->a() * (1 + this->e());
this->r_perigee = this->a() * (1 - this->e());
return;
}
/*! This method only changes the outputs dependent on true anomaly so that one
* orbit may be queried at various points along the orbit*/
Expand All @@ -192,17 +176,14 @@ void KeplerianOrbit::change_f(){
this->velocity_BP_P = cArray2EigenVector3d(v); //
this->true_anomaly_rate = this->n() * pow(this->a(), 2) * sqrt(1 - pow(this->e(), 2)) / pow(this->r(), 2); //
this->radial_rate = this->r() * this->fDot() * this->e() * sin(this->f()) / (1 + this->e() * cos(this->f())); //
this->eccentric_anomaly = safeAcos((this->e() + cos(this->f()) / (1 + this->e() * cos(this->f())))); //
this->eccentric_anomaly = safeAcos(this->e() + cos(this->f()) / (1 + this->e() * cos(this->f()))); //
this->mean_anomaly = this->E() - this->e() * sin(this->E()); //
this->flight_path_angle = safeAcos(sqrt((1 - pow(this->e(), 2)) / (1 - pow(this->e(), 2)*pow(cos(this->E()), 2)))); //
return;
}

/*! This method sets the planet being orbited */
void KeplerianOrbit::set_planet(GravBodyData *plt){
this->planet = plt;
this->mu = plt->mu;
return;
/*! This method sets the gravitational constants of the body being orbited */
void KeplerianOrbit::set_mu(const double mu){
this->mu = mu;
}


Expand Down
93 changes: 46 additions & 47 deletions src/architecture/utilities/keplerianOrbit.h
Original file line number Diff line number Diff line change
Expand Up @@ -20,8 +20,8 @@
#pragma once

#include <Eigen/Dense>
#include "simulation/dynamics/_GeneralModuleFiles/gravityEffector.h"
#include <architecture/utilities/orbitalMotion.h>
#include "architecture/utilities/astroConstants.h"


//! @brief The KeplerianOrbit class represents an elliptical orbit and provides a coherent set of
Expand All @@ -30,37 +30,37 @@
class KeplerianOrbit {
public:
KeplerianOrbit();
KeplerianOrbit(classicElements oe, GravBodyData* planet);
patkenneally marked this conversation as resolved.
Show resolved Hide resolved
KeplerianOrbit(classicElements oe, const double mu);
KeplerianOrbit(const KeplerianOrbit &orig);
~KeplerianOrbit();


Eigen::Vector3d r_BP_P(); //!< body position vector relative to planet
Eigen::Vector3d v_BP_P(); //!< body velocity vector relative to planet
Eigen::Vector3d h_BP_P(); //!< angular momentum of body relative to planet
double M();
double n();
double P();
double f();
double fDot();
double RAAN();
double omega();
double i();
double e();
double a();
double h();
Eigen::Vector3d r_BP_P() const; //!< body position vector relative to planet
Eigen::Vector3d v_BP_P() const; //!< body velocity vector relative to planet
Eigen::Vector3d h_BP_P() const; //!< angular momentum of body relative to planet
double M() const;
double n() const;
double P() const;
double f() const;
double fDot() const;
double RAAN() const;
double omega() const;
double i() const;
double e() const;
double a() const;
double h() const;
double Energy();
double r();
double v();
double r_a();
double r_p();
double fpa();
double E();
double p();
double rDot();
double c3();
double r() const;
double v() const;
double r_a() const;
double r_p() const;
double fpa() const;
double E() const;
double p() const;
double rDot() const;
double c3() const;
classicElements oe();
void set_planet(GravBodyData* plt);
void set_mu(const double mu);
void set_a(double a);
void set_e(double e);
void set_i(double i);
Expand All @@ -69,27 +69,26 @@ class KeplerianOrbit {
void set_f(double f);

private:
GravBodyData* planet;
double mu;
double semi_major_axis;
double eccentricity;
double inclination;
double argument_of_periapsis;
double right_ascension;
double true_anomaly;
double true_anomaly_rate;
double orbital_period;
double orbital_energy;
double v_infinity;
double orbit_radius;
double radial_rate;
double r_apogee;
double r_perigee;
double semi_parameter;
double flight_path_angle;
double eccentric_anomaly;
double mean_motion;
double mean_anomaly;
double mu = MU_EARTH;
double semi_major_axis = 1E5;
double eccentricity = 1E-5;
double inclination{};
double argument_of_periapsis{};
double right_ascension{};
double true_anomaly{};
double true_anomaly_rate{};
double orbital_period{};
double orbital_energy{};
double v_infinity{};
double orbit_radius{};
double radial_rate{};
double r_apogee{};
double r_perigee{};
double semi_parameter{};
double flight_path_angle{};
double eccentric_anomaly{};
double mean_motion{};
double mean_anomaly{};
Eigen::Vector3d orbital_angular_momentum_P;
Eigen::Vector3d position_BP_P;
Eigen::Vector3d velocity_BP_P;
Expand Down
Original file line number Diff line number Diff line change
Expand Up @@ -26,9 +26,7 @@
# Creation Date: Sept 10 2019
patkenneally marked this conversation as resolved.
Show resolved Hide resolved
#

import pytest
import os, inspect
from Basilisk.simulation import gravityEffector
from Basilisk.utilities import orbitalMotion
from Basilisk.architecture import keplerianOrbit
from copy import copy
Expand Down Expand Up @@ -81,9 +79,7 @@ def unitKeplerianOrbit(show_plots=False):

# constructor with arguments
oe = orb.oe()
earth = gravityEffector.GravBodyData()
earth.mu = orbitalMotion.MU_EARTH
orb2 = keplerianOrbit.KeplerianOrbit(oe, earth)
orb2 = keplerianOrbit.KeplerianOrbit(oe, orbitalMotion.MU_EARTH)
assert orb2.r_BP_P() == orb.r_BP_P()
if not orb2.r_BP_P() == orb.r_BP_P():
testFailCount += 1
Expand Down