CARMA C++
Planet.h
Go to the documentation of this file.
1 
11 #ifndef CARMA_SERVICES_PLANET_H
12 #define CARMA_SERVICES_PLANET_H
13 
14 #include "carma/services/Astro.h"
16 #include "carma/services/PlanetTemperature.h"
17 #include "carma/services/Types.h"
18 #include "carma/util/Time.h"
19 #include <memory>
20 #include <string>
21 #include <vector>
22 
23 namespace carma {
24  namespace services {
25 
26  class Angle;
27  class Distance;
28  class Frequency;
29  class FluxDensity;
30  class Length;
31  class Location;
32  class Temperature;
33  class Velocity;
34 
45  class Planet {
46  public:
47 
54  explicit Planet( const std::string& planetName );
55 
57  virtual ~Planet();
58 
59  const std::string getName() const;
60 
69  void setMJD( double mjd = carma::util::Time::MJD() );
70 
75  void setLocation( const Location& location );
76 
82  const Distance earthDistance() const;
83 
87  const Distance avgSunDistance() const;
88 
99  const Angle majorAxis() const;
100 
106  const Angle minorAxis() const;
107 
121  const Angle axisAngle() const;
122 
130  const Angle tiltAngle() const;
131 
143  const Temperature brightnessTemperature(const Frequency& frequency);
144 
150  const FluxDensity flux(const Frequency& frequency);
151 
159  FluxDensity uvFlux(const Frequency& frequency, const Length& u, const Length& v);
160 
167  std::vector<FluxDensity*> uvFluxes(const Frequency& frequency,
168  const std::vector<Length*>& u,
169  const std::vector<Length*>& v);
170 
176 
177 
178  const Table getTable() const {return planetTemperature_->getTable();}
179 
180  private:
183 
184  // the planet name. store here b/c I don't want to
185  // call Ephemeris::getSource().getName() in a catch block.
186  std::string name_;
187 
188  std::auto_ptr<carma::services::PlanetTemperature> planetTemperature_;
189 
191  // This is mutable so that its internal data may change even in
192  // const methods of Planet. I.e. we want to call non-const methods
193  // of Ephemeris inside const methods of Planet.
194  mutable Ephemeris ephemeris_;
195 
196 
197  const Temperature powerLawTb(const Frequency& frequency) const;
198 
208  const Angle planetAngle(bool positionAngle) const;
209 
210  };
211 
212  }
213 }
214 
215 #endif //CARMA_SERVICES_PLANET_H
const Distance avgSunDistance() const
Common time functions.
const Velocity velocity(velocityFrameType frameType)
Astronomical Constants used across CARMA.
void setMJD(double mjd=carma::util::Time::MJD())
Set the Modified Julian Date to use in calculations.
The FluxDensity class is used to represent a flux density in any units.
Definition: FluxDensity.h:20
void setLocation(const Location &location)
Set the location to use in calculations.
const Angle tiltAngle() const
const Angle axisAngle() const
The Distance class can represent an distance in any units.
Definition: Distance.h:43
unsigned int frameType
Half second frames since Jan 1, 2000.
Definition: types.h:29
Various type definitions for services classes.
Planet(const std::string &planetName)
Construct a Planet given the name.
The Velocity class can represent an velocity in any units.
Definition: Velocity.h:42
const Temperature brightnessTemperature(const Frequency &frequency)
The Frequency class can represent any frequency in any units.
Definition: Frequency.h:39
Ephemeris wraps the NOVAS library and any other ephemeris related functions into a simple to use clas...
Definition: Ephemeris.h:91
const Distance earthDistance() const
The Temperature class represents a temperature in any unit.
Definition: Temperature.h:39
const Angle minorAxis() const
A struct to hold simple planetary data.
Definition: Astro.h:55
std::vector< FluxDensity * > uvFluxes(const Frequency &frequency, const std::vector< Length * > &u, const std::vector< Length * > &v)
Simple ASCII Table format, fully memory based.
Definition: Table.h:126
Location specifies a location (observatory if you wish) on planet earth, as longitude, latitude, and altitude above sea-level.
Definition: Location.h:32
const FluxDensity flux(const Frequency &frequency)
The Angle class can represent any angle in any units.
Definition: Angle.h:38
const Angle majorAxis() const
Proposed planet class.
Definition: Planet.h:45
The Length class can represent a length in any units.
Definition: Length.h:36
virtual ~Planet()
Destructor.
enum carma::services::velocityFrameEnum velocityFrameType
The Velocity Frame.
Definition: Types.h:25
FluxDensity uvFlux(const Frequency &frequency, const Length &u, const Length &v)
Declaration of carma::services::Ephemeris.
static double MJD()
Get current MJD.