SatAndLight  2.2.2-hubble
Simulation toolkit for space telescopes
SourceRead.h
Go to the documentation of this file.
1 #ifndef __SourceRead__
2 #define __SourceRead__
3 
10 #include "CUtils.h"
11 #include "SpecTime.h"
12 #include <TChain.h>
13 #include <TH2.h>
14 #include <TGraph.h>
15 #include <TMath.h>
16 
23 class SourceRead: public TChain{
24 
25  public:
26 
32 
38  SourceRead(const string aFilePattern, const string aName);
39 
41 
44  virtual ~SourceRead(void);
50 
53  inline TH2D* GetFlux(void){ return src_sp->GetFlux(); };
54 
56 
59  inline TH2D* GetRate(void){ return src_sp->GetRate(); };
60 
62 
65  inline int GetType(void){ return src_type; };
66 
68 
71  inline TGraph* GetThetap(void){ return src_thetap; };
72 
74 
78  inline double GetThetap(const ULong64_t aTime){ return src_thetap->Eval(aTime); };
79 
81 
84  inline TGraph* GetPhip(void){ return src_phip; };
85 
87 
91  inline double GetPhip(const ULong64_t aTime){ return src_phip->Eval(aTime); };
92 
94 
98  inline double GetRightAscensionDeg(const ULong64_t aTime){ return src_phip->Eval(aTime)*180.0/TMath::Pi(); };
99 
101 
105  inline double GetDeclinationDeg(const ULong64_t aTime){ return 90.0-src_thetap->Eval(aTime)*180.0/TMath::Pi(); };
106 
108 
111  inline TGraph* GetCosConeAngle(void){ return src_cosconeangle; };
112 
114 
118  inline double GetCosConeAngle(const ULong64_t aTime){ return src_cosconeangle->Eval(aTime); };
119 
121 
124  inline string GetName(void){ return src_name; };
125 
127 
131  inline ULong64_t GetAstroTime(void){ return astro_time; };
132 
134 
138  inline double GetAstroEnergy(void){ return astro_energy; };
139 
141 
145  inline double GetAstroThetap(void){ return astro_thetap; };
146 
148 
152  inline double GetAstroPhip(void){ return astro_phip; };
153 
155 
159  inline double GetAstroRightAscensionDeg(void){ return astro_phip*180.0/TMath::Pi(); };
160 
162 
166  inline double GetAstroDeclinationDeg(void){ return 90.0-astro_thetap*180.0/TMath::Pi(); };
167 
169 
172  inline double GetTimeMin(void){ return src_sp->GetTimeMin(); };
173 
175 
178  inline double GetTimeMax(void){ return src_sp->GetTimeMax(); };
179 
181 
187  inline TH1D* GetLightCurve(const double aEnergyMin=0.0, const double aEnergyMax=-1.0){
188  return src_sp->GetLightCurve(aEnergyMin, aEnergyMax);
189  }
190 
192 
198  inline TH1D* GetSpectrum(const double aTimeMin=0.0, const double aTimeMax=-1.0){
199  return src_sp->GetSpectrum(aTimeMin, aTimeMax);
200  }
201 
202  private:
203 
204  // Astroparticles
205  ULong64_t astro_time;
206  double astro_thetap;
207  double astro_phip;
208  double astro_energy;
209 
210  // Source
211  string src_name;
212  int src_type;
213  TGraph *src_thetap;
214  TGraph *src_phip;
217 
218  ClassDef(SourceRead,0)
219 };
220 
221 #endif
222 
223 
This module offers C++ utility functions.
This module is used to describe spectra and light curves of astrophysical sources.
Read a source file sequence.
Definition: SourceRead.h:23
double GetTimeMax(void)
Returns the source time end [ms].
Definition: SourceRead.h:178
string GetName(void)
Returns the source name.
Definition: SourceRead.h:124
TGraph * src_cosconeangle
source cone angle [rad].
Definition: SourceRead.h:215
double GetCosConeAngle(const ULong64_t aTime)
Returns the source cone angle cosine at a given time.
Definition: SourceRead.h:118
virtual ~SourceRead(void)
SourceRead class destructor.
Definition: SourceRead.cc:72
SpecTime * src_sp
source SpecTime object (read-only).
Definition: SourceRead.h:216
TH1D * GetLightCurve(const double aEnergyMin=0.0, const double aEnergyMax=-1.0)
Returns the source light curve between two energies.
Definition: SourceRead.h:187
double GetPhip(const ULong64_t aTime)
Returns the source position angle [rad]$ at a given time.
Definition: SourceRead.h:91
TGraph * GetThetap(void)
Returns the source position angle [rad] as a function of time.
Definition: SourceRead.h:71
double GetAstroThetap(void)
Returns the astro-particle [rad].
Definition: SourceRead.h:145
double GetThetap(const ULong64_t aTime)
Returns the source position angle [rad] at a given time.
Definition: SourceRead.h:78
double GetRightAscensionDeg(const ULong64_t aTime)
Returns the source right ascension in [deg] at a given time.
Definition: SourceRead.h:98
double astro_phip
astroparticle phi prime.
Definition: SourceRead.h:207
double GetAstroDeclinationDeg(void)
Returns the astro-particle declination in [deg].
Definition: SourceRead.h:166
SourceRead(const string aFilePattern, const string aName)
SourceRead class constuctor.
int GetType(void)
Returns the source type (see Source::sourcetype).
Definition: SourceRead.h:65
double astro_thetap
astroparticle .
Definition: SourceRead.h:206
ULong64_t astro_time
astroparticle time [ms].
Definition: SourceRead.h:205
TH2D * GetRate(void)
Returns the source rate.
Definition: SourceRead.h:59
double astro_energy
astroparticle energy [keV].
Definition: SourceRead.h:208
double GetDeclinationDeg(const ULong64_t aTime)
Returns the source declination in [deg] at a given time.
Definition: SourceRead.h:105
double GetAstroRightAscensionDeg(void)
Returns the astro-particle right-ascension in [deg].
Definition: SourceRead.h:159
int src_type
source type, see Source::sourcetype.
Definition: SourceRead.h:212
TGraph * GetCosConeAngle(void)
Returns the source cone angle cosine as a function of time.
Definition: SourceRead.h:111
string src_name
source name.
Definition: SourceRead.h:211
ULong64_t GetAstroTime(void)
Returns the astro-particle time [ms].
Definition: SourceRead.h:131
double GetAstroPhip(void)
Returns the astro-particle [rad].
Definition: SourceRead.h:152
TGraph * src_thetap
source .
Definition: SourceRead.h:213
double GetTimeMin(void)
Returns the source time start [ms].
Definition: SourceRead.h:172
double GetAstroEnergy(void)
Returns the astro-particle energy [keV].
Definition: SourceRead.h:138
TH1D * GetSpectrum(const double aTimeMin=0.0, const double aTimeMax=-1.0)
Returns the source energy spectrum between two times.
Definition: SourceRead.h:198
TGraph * GetPhip(void)
Returns the source position angle [rad]$ as a function of time.
Definition: SourceRead.h:84
TH2D * GetFlux(void)
Returns the source flux.
Definition: SourceRead.h:53
TGraph * src_phip
source .
Definition: SourceRead.h:214
Astrophysical source spectra and light curves.
Definition: SpecTime.h:32
TH1D * GetLightCurve(const double aEnergyMin=0.0, const double aEnergyMax=-1.0)
Returns the source light curve between two energies.
Definition: SpecTime.cc:268
double GetTimeMax(void)
Returns the maximum time defining the SpecTime object.
Definition: SpecTime.h:151
TH2D * GetFlux(void)
Returns a pointer to the flux TH2D object.
Definition: SpecTime.h:387
TH2D * GetRate(void)
Returns a pointer to the flux TH2D object.
Definition: SpecTime.h:394
double GetTimeMin(void)
Returns the minimum time defining the SpecTime object.
Definition: SpecTime.h:157
TH1D * GetSpectrum(const double aTimeMin=0.0, const double aTimeMax=-1.0, const bool aObserved=false)
Returns the source energy spectrum between two times.
Definition: SpecTime.cc:304