balloon.hh
1 #ifndef BALLOON_H
2 #define BALLOON_H
3 
4 
6 //class Balloon:
8 
9 #include <iostream>
10 
11 class Position;
12 class Vector;
13 class TChain;
14 class TTree;
15 class TFile;
16 class TTreeIndex;
17 class IceModel;
18 class Anita;
19 class Settings;
20 class Interaction;
21 class Ray;
22 class TH1F;
23 class TGraph;
24 
25 using std::string;
26 using std::cout;
27 
28 
30 class Balloon {
31 
32 private:
33  // string anitaliteflight; // the gps path of the anita-lite flight
34  //string anitaflight;// gps path of anita flight
35 
36 public:
37  Balloon();
38 
39 
40  // GPS positions of Anita-lite balloon flight
41  int igps;
42  int ibnposition;
43  double dtryingposition;
44  static const int MAX_POSITIONS=50;
45  Vector slacpositions[MAX_POSITIONS];
46  string sslacpositions[MAX_POSITIONS];
47  int islacposition;
48  TChain *flightdatachain;
49  TTreeIndex *tindex;
50  unsigned int realTime_flightdata_temp;
51  unsigned int realTime_flightdata;
52  float flatitude,flongitude,faltitude,fheading,froll, fpitch;
53  double latitude,longitude,altitude,heading,roll,pitch;
54  double MINALTITUDE;
57  int WHICHPATH;
59  double BN_ALTITUDE;
60  unsigned short surfTrigBandMask[9][2];
61  float powerthresh[9][32];
62  float meanp[9][32];
63  double altitude_bn;
64  double theta_bn;
65  double phi_bn;
67  double horizcoord_bn;
68  double vertcoord_bn;
70  Vector x_axis_rot;
71  Vector y_axis_rot;
72  Vector z_axis_rot;
78  double MAXHORIZON;
79  double phi_spin;
80  int NPOINTS;
82  int NPOINTS_MAX;
83  double latitude_bn_anitalite[100000];
84  double longitude_bn_anitalite[100000];
85  double altitude_bn_anitalite[100000];
86  double heading_bn_anitalite[100000];
87  double realtime_bn_anitalite[100000];
88  double BN_LONGITUDE;
89  double BN_LATITUDE;
90  double min_time;
91  double max_time;
92 
93 
95 
105  void setObservationLocation(Interaction *interaction1,int inu,IceModel *antarctic,Settings *settings1);
106 
108 
117  void GetBoresights(Settings *settings1,Anita *anita1,Position r_boresights[Anita::NLAYERS_MAX][Anita::NPHI_MAX]);
118 
120 
132  void PickDownwardInteractionPoint(Interaction *interaction1,Anita *anita1,Settings *settings1,IceModel *antarctica1,
133  Ray *ray1, int &beyondhorizon, double len_int_km2, Vector * force_dir = 0);
134 
136 
142  void InitializeBalloon();
143 
145 
151  void ReadAnitaliteFlight();
152 
154 
161  void CenterPayload(double& hitangle_e);
162 
163 
165 
175  void PickBalloonPosition(Vector straightup,IceModel *antarctica,Settings *settings1, Anita *anita1);
176 
178 
189  void PickBalloonPosition(IceModel *antarctica1,Settings *settings1,int inu,Anita *anita1, double randomNumber);
190 
192 
198  int Getibnposition();
199 
201 
208  double GetBalloonSpin(double heading);
209 
210 
212 
225  void GetAntennaOrientation(Settings *settings1,
226  Anita *anita1,
227  int ilayer, int ifold,
228  Vector& n_eplane,
229  Vector& n_hplane,
230  Vector& n_normal);
231 
233 
246  void GetEcompHcompEvector(Settings *settings1,
247  Vector n_eplane,
248  Vector n_hplane,
249  const Vector n_pol,
250  double& e_component,
251  double& h_component,
252  double& n_component);
253 
255 
268  void GetEcompHcompkvector(Vector n_eplane,
269  Vector n_hplane,
270  Vector n_normal,
271  const Vector n_exit2bn,
272  double& e_component_kvector,
273  double& h_component_kvector,
274  double& n_component_kvector);
275 
277 
288  void GetHitAngles(
289  double e_component_kvector,
290  double h_component_kvector,
291  double n_component_kvector,
292  double& hitangle_e,
293  double& hitangle_h);
294 
296 
303  void SetDefaultBalloonPosition(IceModel *antarctica1);
304 
306 
314  void setr_bn(double latitude,double longitude);
315 
317 
324  void AdjustSlacBalloonPosition(int inu);
325 
327 
334  void GetSlacPositions(Anita *anita1);
335 
337 
345  void GetBoresights(Settings *settings1,Anita *anita1);
346 
348 
356  void calculate_antenna_positions(Settings *settings1,Anita *anita1);
357 
359 
366  Vector RotatePayload(Vector ant_pos);
367 
369 
376  Vector unRotatePayload(Vector ant_pos);
377 
378 
379 
380 }; //class Balloon
381 
382 
383 #endif
384 
385 
unsigned int realTime_flightdata_temp
realtime from the flight data file
Definition: balloon.hh:50
double realtime_bn_anitalite[100000]
same for real life time
Definition: balloon.hh:87
int WHICHPATH
0=fixed balloon position,1=randomized,2=ANITA-lite GPS data,3=banana plot
Definition: balloon.hh:57
unsigned short surfTrigBandMask[9][2]
Ryan&#39;s 16 bit masks for 9 surfs. 2x16 bit masks gives 32 channels per surf.
Definition: balloon.hh:60
int igps
which balloon position do we use out of the 25000 anitalite GPS positions.
Definition: balloon.hh:41
void ReadAnitaliteFlight()
This function reads in the ANITA LITE flight.
Definition: balloon.cc:112
Position r_bn
position of balloon
Definition: balloon.hh:66
static const int MAX_POSITIONS
for the slac beam test
Definition: balloon.hh:44
double phi_bn
theta,phi of balloon wrt south pole
Definition: balloon.hh:65
Vector RotatePayload(Vector ant_pos)
This function rotates the payload.
Definition: balloon.cc:1164
double heading_bn_anitalite[100000]
same for heading of the balloon
Definition: balloon.hh:86
double MAXHORIZON
pick the interaction within this distance from the balloon so that it is within the horizon ...
Definition: balloon.hh:78
int Getibnposition()
This function gets ith balloon position.
Definition: balloon.cc:303
void PickBalloonPosition(Vector straightup, IceModel *antarctica, Settings *settings1, Anita *anita1)
This function picks the balloon position.
Definition: balloon.cc:317
void AdjustSlacBalloonPosition(int inu)
This function adjusts the slac balloon position.
Definition: balloon.cc:627
void CenterPayload(double &hitangle_e)
This function centers the payload.
Definition: balloon.cc:639
double horizcoord_bn
x component of balloon position
Definition: balloon.hh:67
int REDUCEBALLOONPOSITIONS
only take every 100th entry in the flight data file
Definition: balloon.hh:56
double altitude_bn_anitalite[100000]
same for altitude
Definition: balloon.hh:85
int NPOINTS_MIN
min and max index for gps positions we want to include in the simulation (to exclude launch and fall)...
Definition: balloon.hh:81
double BN_LATITUDE
balloon latitude for fixed balloon location
Definition: balloon.hh:89
void PickDownwardInteractionPoint(Interaction *interaction1, Anita *anita1, Settings *settings1, IceModel *antarctica1, Ray *ray1, int &beyondhorizon, double len_int_km2, Vector *force_dir=0)
This function picks downward interaction point.
Definition: balloon.cc:762
Position r_bn_shadow
position of the balloon projected on earth surface - point just below balloon at surface of the earth...
Definition: balloon.hh:77
double MINALTITUDE
minimum altitude balloon needs to be before we consider it a good event to read from the flight data ...
Definition: balloon.hh:54
Reads in and stores input settings for the run.
Definition: Settings.h:37
void setr_bn(double latitude, double longitude)
This function sets r of the balloon.
Definition: balloon.cc:742
int NPOINTS
number of GPS positions we&#39;re picking from.
Definition: balloon.hh:80
double BN_ALTITUDE
pick balloon altitude
Definition: balloon.hh:59
double BN_LONGITUDE
balloon longitude for fixed balloon location
Definition: balloon.hh:88
Vector n_bn
normalized r_bn
Definition: balloon.hh:73
float powerthresh[9][32]
power threshold in Watts
Definition: balloon.hh:61
void InitializeBalloon()
This function initializes the balloon or the specific flight.
Definition: balloon.cc:188
This class is a 3-vector that represents a position on the Earth&#39;s surface.
Definition: position.hh:26
int RANDOMIZE_BN_ORIENTATION
0=fixed balloon orientation,1=randomized
Definition: balloon.hh:58
Stores everything about a particular neutrino interaction. Interaction.
Definition: Primaries.h:135
void GetAntennaOrientation(Settings *settings1, Anita *anita1, int ilayer, int ifold, Vector &n_eplane, Vector &n_hplane, Vector &n_normal)
This function gets the antenna orientation.
Definition: balloon.cc:649
Vector n_north
north, as seen from the balloon position
Definition: balloon.hh:75
double longitude_bn_anitalite[100000]
same for longitude
Definition: balloon.hh:84
double GetBalloonSpin(double heading)
This function gets the spin of the balloon whatever that means.
Definition: balloon.cc:287
double phi_spin
orientation of the balloon
Definition: balloon.hh:79
void SetDefaultBalloonPosition(IceModel *antarctica1)
This function sets the default balloon position.
Definition: balloon.cc:75
Contains everything about positions within payload and signals it sees for each event, in both the trigger and signal paths.
Definition: anita.hh:32
void calculate_antenna_positions(Settings *settings1, Anita *anita1)
This function calculates antenna positions.
Definition: balloon.cc:1127
Vector unRotatePayload(Vector ant_pos)
This function UN-rotates the payload.
Definition: balloon.cc:1221
double surface_under_balloon
distance between center of the earth and the surface of earth under balloon
Definition: balloon.hh:76
static const int NLAYERS_MAX
max number of layers (in smex design, it&#39;s 4)
Definition: anita.hh:59
double vertcoord_bn
y component of balloon position
Definition: balloon.hh:68
double latitude_bn_anitalite[100000]
latitude at times along flightpath, equally distributed among gps data. This is filled with anita or ...
Definition: balloon.hh:83
void setObservationLocation(Interaction *interaction1, int inu, IceModel *antarctic, Settings *settings1)
This function sets the observation location.
Definition: balloon.cc:46
void GetEcompHcompEvector(Settings *settings1, Vector n_eplane, Vector n_hplane, const Vector n_pol, double &e_component, double &h_component, double &n_component)
This function gets the e-component, h-component and e-vector.
Definition: balloon.cc:704
Vector n_east
east, as seen from the balloon position
Definition: balloon.hh:74
Handles everything related to balloon positions, payload orientation over the course of a flight...
Definition: balloon.hh:30
This class represents a three-vector. Operators are overloaded to provide for the familiar operations...
Definition: vector.hh:27
static const int NPHI_MAX
max number of antennas around in phi (in smex, 16)
Definition: anita.hh:61
void GetHitAngles(double e_component_kvector, double h_component_kvector, double n_component_kvector, double &hitangle_e, double &hitangle_h)
This function gets the hit angles.
Definition: balloon.cc:726
int igps_previous
which entry from the flight data file the previous event was so we can just take the next one...
Definition: balloon.hh:55
float meanp[9][32]
mean power in Watts
Definition: balloon.hh:62
double dtryingposition
weighting factor: how many equivalent tries each neutrino counts for after having reduced possible in...
Definition: balloon.hh:43
Ray tracing.
Definition: ray.hh:20
void GetEcompHcompkvector(Vector n_eplane, Vector n_hplane, Vector n_normal, const Vector n_exit2bn, double &e_component_kvector, double &h_component_kvector, double &n_component_kvector)
This function gets the e-component, h-component and k-vector.
Definition: balloon.cc:690
void GetBoresights(Settings *settings1, Anita *anita1, Position r_boresights[Anita::NLAYERS_MAX][Anita::NPHI_MAX])
This function gets the boresights.
Position r_boresights[Anita::NLAYERS_MAX][Anita::NPHI_MAX]
position of antenna boresights
Definition: balloon.hh:69
Ice thicknesses and water depth.
Definition: icemodel.hh:103
unsigned int realTime_flightdata
realtime from the flight data file
Definition: balloon.hh:51
void GetSlacPositions(Anita *anita1)
This function gets the slac balloon positions.
Definition: balloon.cc:917