Skip to content

Commit

Permalink
Add root macro for em profile
Browse files Browse the repository at this point in the history
Add root macro that calculates the longitudinal profile of em showers.
  • Loading branch information
lopezzot committed Jan 24, 2024
1 parent 00bef0c commit 8b7276a
Showing 1 changed file with 75 additions and 0 deletions.
75 changes: 75 additions & 0 deletions analysis/emprofile.C
Original file line number Diff line number Diff line change
@@ -0,0 +1,75 @@
// Simple macro card to reconstruct longitdunal profiles of en-showers

void DoAnalysis();

void emprofile()
{
DoAnalysis();
}

void DoAnalysis()
{
const string filename = "50kFTFP_BERT_EMZ_HGCALTBout_Run0.root";
cout << "Analysis of " << filename << endl;
TFile* file = TFile::Open(filename.c_str(), "READ");
TTree* tree = (TTree*)file->Get("HGCALTBout");

double CEETot{0.};
tree->SetBranchAddress("CEETot", &CEETot);
double CHETot{0.};
tree->SetBranchAddress("CHETot", &CHETot);
double AHCALTot{0.};
tree->SetBranchAddress("AHCALTot", &AHCALTot);
double HGCALTot{0.};
tree->SetBranchAddress("HGCALTot", &HGCALTot);
int IntLayer{0};
tree->SetBranchAddress("IntLayer", &IntLayer);
vector<double>* CEESignals = NULL;
tree->SetBranchAddress("CEESignals", &CEESignals);
vector<double>* CHESignals = NULL;
tree->SetBranchAddress("CHESignals", &CHESignals);
vector<double>* AHCALSignals = NULL;
tree->SetBranchAddress("AHCALSignals", &AHCALSignals);

const int layersNo = 28;
std::array<double, layersNo/2> emprofile{0.};
std::array<double, layersNo> fullemprofile{0.};

for (std::size_t evtNo = 0; evtNo < tree->GetEntries(); evtNo++) {
tree->GetEntry(evtNo);

for(std::size_t i = 0; i < layersNo; i = i+2){
emprofile.at(i/2) += CEESignals->at(i)+CEESignals->at(i+1);
}
for(std::size_t i = 0; i < layersNo; i++){
fullemprofile.at(i) += CEESignals->at(i);
}
}

for(const auto& a : emprofile) cout<<a<<endl;
cout<<"-----------------------"<<endl;
for(auto& a : emprofile) { a = a / tree->GetEntries(); }
for(const auto& a : emprofile) cout<<a<<endl;

for(auto& a : fullemprofile) { a = a / tree->GetEntries(); }

std::array<double, layersNo/2> layers{};
for (std::size_t i=0; i<layersNo/2; i++) { layers.at(i) = i; }
std::array<double, layersNo> fulllayers{};
for (std::size_t i=0; i<layersNo; i++) { fulllayers.at(i) = i; }

TGraph grEmProfile(layersNo/2, layers.data(), emprofile.data());
grEmProfile.SetTitle("EmProfile");
grEmProfile.SetName("EmProfile");

TGraph grFullEmProfile(layersNo, fulllayers.data(), fullemprofile.data());
grFullEmProfile.SetTitle("FullEmProfile");
grFullEmProfile.SetName("FullEmProfile");

TFile* outputfile(TFile::Open("HGCALTBemprofile.root", "RECREATE"));
outputfile->cd();
grEmProfile.Write();
grFullEmProfile.Write();
outputfile->Close();

}

0 comments on commit 8b7276a

Please sign in to comment.