SatAndLight  2.2.2-hubble
Simulation toolkit for space telescopes
SpecTime.h
Go to the documentation of this file.
1 #ifndef __SpecTime__
2 #define __SpecTime__
3 
9 #include "CUtils.h"
10 #include <TMath.h>
11 #include <TF1.h>
12 #include <TH2.h>
13 #include <TGraphAsymmErrors.h>
14 #include <TRandom3.h>
15 
16 using namespace std;
17 
32 class SpecTime{
33 
34  public:
35 
41 
50  SpecTime(const int aNt, const double *aTimeBins,
51  const int aNe, const double *aEnergyBins);
52 
54 
66  SpecTime(const int aNt, const double aTimeMin, const double aTimeMax, const bool aUseTimeLog,
67  const int aNe, const double aEnergyMin, const double aEnergyMax, const bool aUseEnergyLog);
68 
70 
78  SpecTime(TH2D* aFlux, TH2D* aRate=NULL);
79 
81 
84  virtual ~SpecTime(void);
90 
96  inline void SetPhotoElectricCrossSection(TGraph *aGraph){
97  delete photoelectric_xsec;
98  photoelectric_xsec = (TGraph*)aGraph->Clone("sp_xsec");
99  photoelectric_xsec->SetBit(TGraph::kIsSortedX);
100  };
101 
103 
113  void SetFlux(TH2D *aFlux, TGraph *aSurface=NULL);
114 
116 
129  void MakeRate(TGraph *aSurface);
130 
132 
145  int GetParticleN(const ULong64_t aTime, const UInt_t aDuration, const double aCosTheta=1.0);
146 
148 
151  inline double GetTimeMax(void){ return flux->GetXaxis()->GetBinUpEdge(flux->GetNbinsX()); };
152 
154 
157  inline double GetTimeMin(void){ return flux->GetXaxis()->GetBinLowEdge(1); };
158 
160 
172  TH1D* GetLightCurve(const double aEnergyMin=0.0, const double aEnergyMax=-1.0);
173 
175 
193  TH1D* GetSpectrum(const double aTimeMin=0.0, const double aTimeMax=-1.0, const bool aObserved=false);
194 
196 
208  void MakeUniform(const double aFlux, TGraph *aSurface=NULL);
209 
211 
229  void MakeBandLimited(const double aFlux,
230  const double aEnergyMin, const double aEnergyMax,
231  TGraph *aSurface=NULL);
232 
234 
254  void MakePowerLawTime(const double aFlux, const double aTimeDecay,
255  const double aEnergyMin, const double aEnergyMax,
256  TGraph *aSurface=NULL);
257 
259 
276  void MakePowerLawEnergy(const double aFlux, const double aEnergyDecay,
277  TGraph *aSurface=NULL);
278 
280 
303  void MakeGRBModel1(const double aFlux, const double aEnergy0,
304  const double aAlpha, const double aBeta,
305  TGraph *aSurface=NULL);
306 
308 
366  void MakeGRBData(TH1D *aLightCurveExp,
367  const double aE1, const double aE2,
368  const double aPhotonIndex,
369  const double aNint, const double aNgal,
370  const double aZ, TGraph *aSurfaceExp,
371  TGraph *aSurface=NULL);
372 
374 
380  double DrawEnergy(const double aTime);
381 
383 
387  inline TH2D* GetFlux(void){ return flux; };
388 
390 
394  inline TH2D* GetRate(void){ return rate; };
395 
396  protected:
397 
398  TH2D *flux;
399  TH2D *rate;
400  TRandom3 *randgen;
401 
402  private:
404  double **rate_cum;
405  void MakeCum(void);
406 
407  void ConstructHistos(const int aNt, const double *aTimeBins,
408  const int aNe, const double *aEnergyBins);
409 
410  ClassDef(SpecTime,0)
411 };
412 
413 #endif
414 
415 
This module offers C++ utility functions.
Astrophysical source spectra and light curves.
Definition: SpecTime.h:32
double GetTimeMax(void)
Returns the maximum time defining the SpecTime object.
Definition: SpecTime.h:151
TH2D * rate
source rate .
Definition: SpecTime.h:399
void ConstructHistos(const int aNt, const double *aTimeBins, const int aNe, const double *aEnergyBins)
for constructors.
TRandom3 * randgen
random generator.
Definition: SpecTime.h:400
TH2D * GetFlux(void)
Returns a pointer to the flux TH2D object.
Definition: SpecTime.h:387
double ** rate_cum
cumulative rate (over ).
Definition: SpecTime.h:404
TH2D * GetRate(void)
Returns a pointer to the flux TH2D object.
Definition: SpecTime.h:394
void SetPhotoElectricCrossSection(TGraph *aGraph)
Sets the photo-electric cross-section.
Definition: SpecTime.h:96
double GetTimeMin(void)
Returns the minimum time defining the SpecTime object.
Definition: SpecTime.h:157
TH2D * flux
source flux .
Definition: SpecTime.h:394
TGraph * photoelectric_xsec
photo-electric cross-section.
Definition: SpecTime.h:403