13 #include <TGraphAsymmErrors.h>
50 SpecTime(
const int aNt,
const double *aTimeBins,
51 const int aNe,
const double *aEnergyBins);
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);
78 SpecTime(TH2D* aFlux, TH2D* aRate=NULL);
97 delete photoelectric_xsec;
98 photoelectric_xsec = (TGraph*)aGraph->Clone(
"sp_xsec");
99 photoelectric_xsec->SetBit(TGraph::kIsSortedX);
113 void SetFlux(TH2D *aFlux, TGraph *aSurface=NULL);
129 void MakeRate(TGraph *aSurface);
145 int GetParticleN(
const ULong64_t aTime,
const UInt_t aDuration,
const double aCosTheta=1.0);
151 inline double GetTimeMax(
void){
return flux->GetXaxis()->GetBinUpEdge(flux->GetNbinsX()); };
157 inline double GetTimeMin(
void){
return flux->GetXaxis()->GetBinLowEdge(1); };
172 TH1D* GetLightCurve(
const double aEnergyMin=0.0,
const double aEnergyMax=-1.0);
193 TH1D* GetSpectrum(
const double aTimeMin=0.0,
const double aTimeMax=-1.0,
const bool aObserved=
false);
208 void MakeUniform(
const double aFlux, TGraph *aSurface=NULL);
229 void MakeBandLimited(
const double aFlux,
230 const double aEnergyMin,
const double aEnergyMax,
231 TGraph *aSurface=NULL);
254 void MakePowerLawTime(
const double aFlux,
const double aTimeDecay,
255 const double aEnergyMin,
const double aEnergyMax,
256 TGraph *aSurface=NULL);
276 void MakePowerLawEnergy(
const double aFlux,
const double aEnergyDecay,
277 TGraph *aSurface=NULL);
303 void MakeGRBModel1(
const double aFlux,
const double aEnergy0,
304 const double aAlpha,
const double aBeta,
305 TGraph *aSurface=NULL);
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);
380 double DrawEnergy(
const double aTime);
408 const int aNe,
const double *aEnergyBins);
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