// $Id: ShowerStep.cc,v 1.11 2005/01/07 22:48:13 thea Exp $
// Author: Alessandro Thea, Sergio Bottai, Dmitry Naumov 2004/12/08
/*****************************************************************************
* ESAF: Euso Simulation and Analysis Framework *
* *
* Id: ShowerStep *
* Package: Shower *
* Coordinator: Sergio.Bottai, Dmitry.Naumov *
* *
*****************************************************************************/
//_____________________________________________________________________________
//
// ShowerStep is the universal container of Atmospheric Air Shower
// produced by Ultra High Energy Cosmic Ray entering the atmosphere.
// ShowerStep is a vector element of the ShowerTrack object which is
// produced by interfaces to specific generators like CORSICA, AIRES, UNISIM, SLAST, etc
// ShowerStep contains relevant information of EAS development at the given
// position in space
//
//
#include "ShowerStep.hh"
#include "ShowerTrack.hh"
#include <TF12.h>
ClassImp(ShowerStep)
//______________________________________________________________________________
ShowerStep::ShowerStep() : fXi(0.),fXf(0.),fXYZi(0.,0.,0.),fXYZf(0.,0.,0.),
fTimei(0.),fTimef(0.),fAgei(0.),fAgef(0.),fNelectrons(0.),fStepID(-1) {
//
// Constructor
//
fNcharged = 0;
fEloss = 0;
fNcherenkov = 0;
fHistoEnergy = NULL;
fHistoLateral = NULL;
fHistoEneAng = NULL;
fHistoRadPhiEle = NULL;
fHistoRadDTimeEle = NULL;
fHistoRadPhiEloss = NULL;
fHistoAngCher = NULL;
fParentTrack = NULL;
}
//______________________________________________________________________________
ShowerStep::ShowerStep(const ShowerStep &s) {
//
// Copy Constructor
//
fXi = s.fXi;
fXf = s.fXf;
fXYZi = s.fXYZi;
fXYZf = s.fXYZf;
fTimei = s.fTimei;
fTimef = s.fTimef;
fAgei = s.fAgei;
fAgef = s.fAgef;
fNelectrons = s.fNelectrons;
fNcharged = s.fNcharged;
fEloss = s.fEloss;
fNcherenkov = s.fNcherenkov;
fParentTrack = s.fParentTrack;
fStepID = fStepID;
fHistoEnergy = s.fHistoEnergy ? (TH1F*)s.fHistoEnergy->Clone() : 0;
fHistoLateral = s.fHistoLateral ? (TH1F*)s.fHistoLateral->Clone() : 0;
fHistoEneAng = s.fHistoEneAng ? (TH2F*)s.fHistoEneAng->Clone() : 0;
fHistoRadPhiEle = s.fHistoRadPhiEle ? (TH2F*)s.fHistoRadPhiEle->Clone() : 0;
fHistoRadDTimeEle = s.fHistoRadDTimeEle ? (TH2F*)s.fHistoRadDTimeEle->Clone() : 0;
fHistoRadPhiEloss = s.fHistoRadPhiEloss ? (TH2F*)s.fHistoRadPhiEloss->Clone() : 0;
fHistoAngCher = s.fHistoAngCher ? (TH1F*)s.fHistoAngCher->Clone() : 0;
}
//______________________________________________________________________________
ShowerStep::~ShowerStep() {
//
// Destructor
//
Clear();
}
//______________________________________________________________________________
void ShowerStep::Clear() {
//
// Clear the step
//
fXi = 0.;
fXf = 0.;
fXYZi.SetXYZ(0.,0.,0.);
fXYZf.SetXYZ(0.,0.,0.);
fTimei = 0.;
fTimef = 0.;
fAgei = 0.;
fAgef = 0.;
fNelectrons = 0;
fStepID =-1;
fNcharged = 0;
fEloss = 0;
fNcherenkov = 0;
fParentTrack = NULL;
SafeDelete(fHistoEnergy);
SafeDelete(fHistoLateral);
SafeDelete(fHistoEneAng);
SafeDelete(fHistoRadPhiEle);
SafeDelete(fHistoRadDTimeEle);
SafeDelete(fHistoRadPhiEloss);
SafeDelete(fHistoAngCher);
}
//______________________________________________________________________________
const TH1F* ShowerStep::GetHistoEnergy() {
//
// This method returns pointer to the histogram with energy distribution of
// electrons at the given step filled by the Shower Generator. If Shower
// Generator does not provide this information the histogram is filled
// automatically according to configured analitical distribution.
//
if (fHistoEnergy)
return fHistoEnergy;
else {
// no parent, no histo
if ( !fParentTrack ) return 0;
// check if the histo for the given age interval already exists. In
// this case use it otherwise create the new one
Float_t age = (fAgei + fAgef)/2;
Int_t Age_index = -1;
Int_t fNumEnergyHistos = fParentTrack->GetNumEnergyHistos();
for (Int_t i(0); i<fNumEnergyHistos; i++) {
if (i*2./fNumEnergyHistos < age && age <= (i+1)*2./fNumEnergyHistos) Age_index = i;
}
// define histo unique name
char name[100], title[100];
sprintf(name,"ShowerStep.HistoEnergy.AgeIndex%d",Age_index);
sprintf(title,"Energy distribution");
Int_t Nx(1000);
Double_t Xbins[Nx+1];
Double_t MaxElectrons(1.e4);
Xbins[0] = 0.;
for (Int_t i(1);i<=Nx;i++) {
if (i<=100)
Xbins[i] = i*0.01; // step size of 10 keV for [0,1] MeV
else if (i>=100 && i<= 500)
Xbins[i] = Xbins[i-1] + 99/400; // step size of 0.25 MeV for [1,100] MeV
else
Xbins[i] = Xbins[i-1] + 1.; // step size of 1 MeV for [100, infty] MeV
}
fHistoEnergy = new TH1F(name,title,Nx,&Xbins[0]);
fHistoEnergy->SetDirectory(0);
TF12 EnergyDistribution("EnergyDistr",fParentTrack->GetEnergyDistribution("default"),
(fAgei + fAgef)/2,"X");
if (fNelectrons<MaxElectrons) {
fHistoEnergy->FillRandom("EnergyDistr",(Long64_t)fNelectrons);
if (fNelectrons) fHistoEnergy->Scale(1/fNelectrons);
}
else {
fHistoEnergy->FillRandom("EnergyDistr",(Long64_t)MaxElectrons);
if (MaxElectrons) fHistoEnergy->Scale(1/MaxElectrons);
}
return fHistoEnergy;
}
}
//______________________________________________________________________________
const TH1F* ShowerStep::GetHistoEnergy() const {
//
// const version of GetHistoEnergy
//
return ((ShowerStep*)this)->GetHistoEnergy();
}
//______________________________________________________________________________
const TH1F* ShowerStep::GetHistoLateral() {
//
// This method returns pointer to the histogram with lateral distribution
// of electrons at the given step filled by the Shower Generator. If Shower
// Generator does not provide this information the histogram is filled
// automatically according to configured analitical distribution.
//
if (fHistoLateral)
return fHistoLateral;
else {
// no parent, no histo
if ( !fParentTrack ) return 0;
// check if the histo for the given age interval already exists. In
// this case use it otherwise create the new one
Float_t age = (fAgei + fAgef)/2;
Int_t Age_index = -1;
Int_t fNumLateralHistos = fParentTrack->GetNumLateralHistos();
for (Int_t i(0); i<fNumLateralHistos; i++) {
if (i*2./fNumLateralHistos < age && age <= (i+1)*2./fNumLateralHistos) Age_index = i;
}
// define histo unique name
char name[100], title[100];
sprintf(name,"ShowerStep.HistoLateral.AgeIndex%d",Age_index);
sprintf(title,"Lateral distribution");
Int_t Nx(1000);
Double_t MaxElectrons(1.e4);
fHistoLateral = new TH1F(name,title,Nx,0,100);
fHistoLateral->SetDirectory(0);
TF12 LateralDistribution("LateralDistr",fParentTrack->GetLateralDistribution("default"),(fAgei + fAgef)/2,"X");
if (fNelectrons<MaxElectrons) {
fHistoLateral->FillRandom("LateralDistr",(Long64_t)fNelectrons);
if (fNelectrons) fHistoLateral->Scale(1/fNelectrons);
}
else {
fHistoLateral->FillRandom("LateralDistr",(Long64_t)MaxElectrons);
if (MaxElectrons) fHistoLateral->Scale(1/MaxElectrons);
}
return fHistoLateral;
}
}
//______________________________________________________________________________
const TH1F* ShowerStep::GetHistoLateral() const {
//
// const version of GetHistoLateral
//
return ((ShowerStep*)this)->GetHistoLateral();
}