Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:10:08

0001 // This file is part of the Acts project.
0002 //
0003 // Copyright (C) 2021 CERN for the benefit of the Acts project
0004 //
0005 // This Source Code Form is subject to the terms of the Mozilla Public
0006 // License, v. 2.0. If a copy of the MPL was not distributed with this
0007 // file, You can obtain one at http://mozilla.org/MPL/2.0/.
0008 
0009 #include <cmath>
0010 #include <iostream>
0011 #include <map>
0012 #include <string>
0013 #include <tuple>
0014 
0015 #include <TCanvas.h>
0016 #include <TDirectory.h>
0017 #include <TFile.h>
0018 #include <TH1F.h>
0019 #include <TProfile.h>
0020 #include <TTree.h>
0021 
0022 struct MaterialHistograms {
0023   TProfile* x0_vs_eta = nullptr;
0024   TProfile* l0_vs_eta = nullptr;
0025 
0026   TProfile* x0_vs_phi = nullptr;
0027   TProfile* l0_vs_phi = nullptr;
0028 
0029   float s_x0 = 0.;
0030   float s_l0 = 0.;
0031 
0032   MaterialHistograms() = default;
0033 
0034   /// Material histogram constructor
0035   ///
0036   /// @param iA the atomic number
0037   MaterialHistograms(const std::string& name, unsigned int iA,
0038                      unsigned int bins, float eta) {
0039     std::string x0NameEta =
0040         (iA == 0) ? name + std::string("_x0_vs_eta_all")
0041                   : name + std::string("_x0_vs_eta_A") + std::to_string(iA);
0042     std::string l0NameEta =
0043         (iA == 0) ? name + std::string("_l0_vs_eta_all")
0044                   : name + std::string("_l0_vs_eta_A") + std::to_string(iA);
0045 
0046     x0_vs_eta = new TProfile(x0NameEta.c_str(), "X_{0} vs. #eta", bins, -eta,
0047                              eta);
0048     l0_vs_eta = new TProfile(l0NameEta.c_str(), "L_{0} vs. #eta", bins, -eta,
0049                              eta);
0050 
0051     std::string x0NamePhi =
0052         (iA == 0) ? name + std::string("_x0_vs_phi_all")
0053                   : name + std::string("_x0_vs_phi_A") + std::to_string(iA);
0054     std::string l0NamePhi =
0055         (iA == 0) ? name + std::string("_l0_vs_phi_all")
0056                   : name + std::string("_l0_vs_phi_A") + std::to_string(iA);
0057 
0058     x0_vs_phi = new TProfile(x0NamePhi.c_str(), "X_{0} vs. #phi", bins, -M_PI,
0059                              M_PI);
0060     l0_vs_phi = new TProfile(l0NamePhi.c_str(), "L_{0} vs. #phi", bins, -M_PI,
0061                              M_PI);
0062   }
0063 
0064   /// This fills the event into the histograms
0065   /// and clears the cache accordingly
0066   ///
0067   /// @param eta the pseudorapidity value
0068   /// @param phi the phi value
0069   ///
0070   void fillAndClear(float eta, float phi) {
0071     x0_vs_eta->Fill(eta, s_x0);
0072     l0_vs_eta->Fill(eta, s_l0);
0073 
0074     x0_vs_phi->Fill(phi, s_x0);
0075     l0_vs_phi->Fill(phi, s_l0);
0076 
0077     s_x0 = 0.;
0078     s_l0 = 0.;
0079   }
0080 
0081   /// Write out the histograms, the TDirectory needs
0082   /// to be set before
0083   ///
0084   /// Histograms with no contribution will not be
0085   /// written to file.
0086   void write() {
0087     if (x0_vs_eta->GetMaximum() > 0.) {
0088       x0_vs_eta->Write();
0089       l0_vs_eta->Write();
0090 
0091       x0_vs_phi->Write();
0092       l0_vs_phi->Write();
0093     }
0094   }
0095 };
0096 
0097 struct Region {
0098   std::string name;
0099   std::vector<std::tuple<float, float, float, float>> boxes;
0100 
0101   bool inside(float r, float z) const {
0102     for (const auto& [minR, maxR, minZ, maxZ] : boxes) {
0103       if (minR <= r && r < maxR && minZ <= z && z < maxZ) {
0104         return true;
0105       }
0106     }
0107     return false;
0108   }
0109 };
0110 
0111 /// Plot the material composition
0112 ///
0113 /// @param inFile the input root file
0114 /// @param treeNAme the input tree name (default: 'trackstates)
0115 /// @param outFile the output root file
0116 /// @param bins the number of bins
0117 /// @param eta the eta range
0118 void materialComposition(const std::string& inFile, const std::string& treeName,
0119                          const std::string& outFile, unsigned int bins,
0120                          float eta, const std::vector<Region>& regions) {
0121   // Open the input file & get the tree
0122   auto inputFile = TFile::Open(inFile.c_str());
0123   auto inputTree = dynamic_cast<TTree*>(inputFile->Get(treeName.c_str()));
0124   if (inputTree == nullptr) {
0125     return;
0126   }
0127 
0128   // Get the different atomic numbers
0129   TCanvas* materialCanvas =
0130       new TCanvas("materialCanvas", "Materials", 100, 100, 620, 400);
0131   materialCanvas->cd();
0132   // Draw all the atomic elements & get the histogram
0133   inputTree->Draw("mat_A>>hA(200,0.5,200.5)");
0134   TH1F* histA = dynamic_cast<TH1F*>(gDirectory->Get("hA"));
0135   histA->Draw();
0136 
0137   auto outputFile = TFile::Open(outFile.c_str(), "recreate");
0138 
0139   float v_eta = 0;
0140   float v_phi = 0;
0141   std::vector<float>* stepLength = new std::vector<float>;
0142   std::vector<float>* stepX0 = new std::vector<float>;
0143   std::vector<float>* stepL0 = new std::vector<float>;
0144   std::vector<float>* stepA = new std::vector<float>;
0145 
0146   std::vector<float>* stepX = new std::vector<float>;
0147   std::vector<float>* stepY = new std::vector<float>;
0148   std::vector<float>* stepZ = new std::vector<float>;
0149 
0150   inputTree->SetBranchAddress("v_eta", &v_eta);
0151   inputTree->SetBranchAddress("v_phi", &v_phi);
0152   inputTree->SetBranchAddress("mat_step_length", &stepLength);
0153 
0154   inputTree->SetBranchAddress("mat_X0", &stepX0);
0155   inputTree->SetBranchAddress("mat_L0", &stepL0);
0156   inputTree->SetBranchAddress("mat_A", &stepA);
0157 
0158   inputTree->SetBranchAddress("mat_x", &stepX);
0159   inputTree->SetBranchAddress("mat_y", &stepY);
0160   inputTree->SetBranchAddress("mat_z", &stepZ);
0161 
0162   // Loop over all entries ---------------
0163   unsigned int entries = inputTree->GetEntries();
0164 
0165 #ifdef BOOST_AVAILABLE
0166   std::cout << "*** Event Loop: " << std::endl;
0167   progress_display event_loop_progress(entries * regions.size());
0168 #endif
0169 
0170   // Loop of the regions
0171   for (auto& region : regions) {
0172     const auto rName = region.name;
0173 
0174     // The material histograms ordered by atomic mass
0175     std::map<unsigned int, MaterialHistograms> mCache;
0176 
0177     // The material histograms for all
0178     mCache[0] = MaterialHistograms(rName, 0, bins, eta);
0179     for (unsigned int ib = 1; ib <= 200; ++ib) {
0180       if (histA->GetBinContent(ib) > 0.) {
0181         mCache[ib] = MaterialHistograms(rName, ib, bins, eta);
0182       }
0183     }
0184 
0185     for (unsigned int ie = 0; ie < entries; ++ie) {
0186       // Get the entry
0187       inputTree->GetEntry(ie);
0188 
0189 #ifdef BOOST_AVAILABLE
0190       ++event_loop_progress;
0191 #endif
0192 
0193       // Accumulate the material per track
0194       std::size_t steps = stepLength->size();
0195       for (unsigned int is = 0; is < steps; ++is) {
0196         float x = stepX->at(is);
0197         float y = stepY->at(is);
0198         float z = stepZ->at(is);
0199         float r = std::hypot(x, y);
0200 
0201         if (!region.inside(r, z)) {
0202           continue;
0203         }
0204 
0205         float step = stepLength->at(is);
0206         float X0 = stepX0->at(is);
0207         float L0 = stepL0->at(is);
0208 
0209         // The integral one
0210         auto& all = mCache[0];
0211         all.s_x0 += step / X0;
0212         all.s_l0 += step / L0;
0213 
0214         unsigned int sA = histA->FindBin(stepA->at(is));
0215         // The current one
0216         auto currentIt = mCache.find(sA);
0217         if (currentIt == mCache.end()) {
0218           throw std::runtime_error{"Unknown atomic number " +
0219                                    std::to_string(sA)};
0220         }
0221         auto& current = currentIt->second;
0222         current.s_x0 += step / X0;
0223         current.s_l0 += step / L0;
0224       }
0225       // Fill the profiles and clear the cache
0226       for (auto& [key, cache] : mCache) {
0227         cache.fillAndClear(v_eta, v_phi);
0228       }
0229     }
0230     // -----------------------------------
0231 
0232     for (auto [key, cache] : mCache) {
0233       cache.write();
0234     }
0235   }
0236 
0237   outputFile->Close();
0238 
0239   delete stepLength;
0240   delete stepX0;
0241   delete stepL0;
0242   delete stepA;
0243 
0244   delete stepX;
0245   delete stepY;
0246   delete stepZ;
0247 
0248   delete materialCanvas;
0249 }