Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2026-05-23 08:10:28

0001 #ifndef V0_DUPLICATE_READER_H
0002 #define V0_DUPLICATE_READER_H
0003 
0004 #include <TTree.h>
0005 
0006 #include <cmath>
0007 #include <set>
0008 #include <stdexcept>
0009 #include <string>
0010 #include <vector>
0011 
0012 class V0DuplicateReader
0013 {
0014  public:
0015   enum class ParticleType
0016   {
0017     K0s,
0018     Lambda
0019   };
0020 
0021   V0DuplicateReader(TTree* tree, ParticleType particleType)
0022   : m_tree(tree)
0023   , m_particleType(particleType)
0024   {
0025     if (!m_tree)
0026     {
0027       throw std::runtime_error("V0DuplicateReader: input tree is null");
0028     }
0029 
0030     setupCommonBranches();
0031     setupParticleBranches();
0032   }
0033 
0034   void enableDeltaBCOCut(long long minDeltaBCO, long long maxDeltaBCO)
0035   {
0036     m_useDeltaBCOCut = true;
0037     m_minDeltaBCO = minDeltaBCO;
0038     m_maxDeltaBCO = maxDeltaBCO;
0039   }
0040 
0041   void disableDeltaBCOCut()
0042   {
0043     m_useDeltaBCOCut = false;
0044   }
0045 
0046   bool loadEntry(Long64_t entry)
0047   {
0048     if (!m_tree) return false;
0049 
0050     m_tree->GetEntry(entry);
0051 
0052     m_deltaBCO = (m_bco >= 0 && m_eventBCO >= 0) ? (m_bco - m_eventBCO) : -1;
0053 
0054     if (m_useDeltaBCOCut)
0055     {
0056       m_passesDeltaBCO = (m_deltaBCO >= m_minDeltaBCO && m_deltaBCO < m_maxDeltaBCO);
0057     }
0058     else
0059     {
0060       m_passesDeltaBCO = true;
0061     }
0062 
0063     if (!m_passesDeltaBCO)
0064     {
0065       m_currentEntryIsUnique = false;
0066       return true;
0067     }
0068 
0069     std::vector<float> candidateVariables = {
0070       m_mass,
0071       m_pt,
0072       m_eta,
0073       m_phi,
0074       m_decayLength,
0075       m_dira,
0076       m_ip,
0077       m_x,
0078       m_y,
0079       m_z,
0080       m_chi2
0081     };
0082 
0083     const std::string eventKey =
0084       std::to_string(m_runNumber) + "_" + std::to_string(m_bco);
0085 
0086     const std::string candidateKey = buildCandidateKey(candidateVariables);
0087     const std::string fullKey = eventKey + "__" + candidateKey;
0088 
0089     const auto [it, inserted] = m_seenCandidateKeys.insert(fullKey);
0090 
0091     if (inserted)
0092     {
0093       ++m_numberOfUniqueCandidates;
0094       m_currentEntryIsUnique = true;
0095     }
0096     else
0097     {
0098       ++m_numberOfDuplicateCandidates;
0099       m_currentEntryIsUnique = false;
0100     }
0101 
0102     return true;
0103   }
0104 
0105   bool isCurrentEntryUnique() const
0106   {
0107     return m_currentEntryIsUnique;
0108   }
0109 
0110   bool passesDeltaBCOCut() const
0111   {
0112     return m_passesDeltaBCO;
0113   }
0114 
0115   long long deltaBCO() const
0116   {
0117     return m_deltaBCO;
0118   }
0119 
0120   void reset()
0121   {
0122     m_seenCandidateKeys.clear();
0123     m_numberOfUniqueCandidates = 0;
0124     m_numberOfDuplicateCandidates = 0;
0125     m_currentEntryIsUnique = false;
0126     m_deltaBCO = -1;
0127     m_passesDeltaBCO = true;
0128   }
0129 
0130   Long64_t entries() const
0131   {
0132     return m_tree ? m_tree->GetEntries() : 0;
0133   }
0134 
0135   int runNumber() const { return m_runNumber; }
0136   Long64_t bco() const { return m_bco; }
0137   Long64_t eventBCO() const { return m_eventBCO; }
0138   int eventNumber() const { return m_eventNumber; }
0139 
0140   float mass() const { return m_mass; }
0141   float pt() const { return m_pt; }
0142   float eta() const { return m_eta; }
0143   float phi() const { return m_phi; }
0144   float decayLength() const { return m_decayLength; }
0145   float dira() const { return m_dira; }
0146   float ip() const { return m_ip; }
0147   float x() const { return m_x; }
0148   float y() const { return m_y; }
0149   float z() const { return m_z; }
0150   float chi2() const { return m_chi2; }
0151   float rapidity() const { return m_rapidity; }
0152 
0153   long long numberOfUniqueCandidates() const
0154   {
0155     return m_numberOfUniqueCandidates;
0156   }
0157 
0158   long long numberOfDuplicateCandidates() const
0159   {
0160     return m_numberOfDuplicateCandidates;
0161   }
0162 
0163  private:
0164   void setupCommonBranches()
0165   {
0166     m_tree->SetBranchAddress("runNumber", &m_runNumber);
0167     m_tree->SetBranchAddress("eventNumber", &m_eventNumber);
0168     m_tree->SetBranchAddress("BCO", &m_bco);
0169     m_tree->SetBranchAddress("event_bco", &m_eventBCO);
0170   }
0171 
0172   void setupParticleBranches()
0173   {
0174     const std::string prefix = (m_particleType == ParticleType::K0s) ? "K_S0" : "Lambda0";
0175 
0176     m_tree->SetBranchAddress((prefix + "_mass").c_str(), &m_mass);
0177     m_tree->SetBranchAddress((prefix + "_pT").c_str(), &m_pt);
0178     m_tree->SetBranchAddress((prefix + "_pseudorapidity").c_str(), &m_eta);
0179     m_tree->SetBranchAddress((prefix + "_phi").c_str(), &m_phi);
0180     m_tree->SetBranchAddress((prefix + "_decayLength").c_str(), &m_decayLength);
0181     m_tree->SetBranchAddress((prefix + "_DIRA").c_str(), &m_dira);
0182     m_tree->SetBranchAddress((prefix + "_IP").c_str(), &m_ip);
0183     m_tree->SetBranchAddress((prefix + "_x").c_str(), &m_x);
0184     m_tree->SetBranchAddress((prefix + "_y").c_str(), &m_y);
0185     m_tree->SetBranchAddress((prefix + "_z").c_str(), &m_z);
0186     m_tree->SetBranchAddress((prefix + "_chi2").c_str(), &m_chi2);
0187     m_tree->SetBranchAddress((prefix + "_rapidity").c_str(), &m_rapidity);
0188   }
0189 
0190   static std::string buildCandidateKey(const std::vector<float>& candidateVariables)
0191   {
0192     std::string key;
0193     key.reserve(candidateVariables.size() * 16);
0194 
0195     for (float value : candidateVariables)
0196     {
0197       const long long roundedValue = llround(value * 1e6);
0198       key += std::to_string(roundedValue);
0199       key += "_";
0200     }
0201 
0202     return key;
0203   }
0204 
0205  private:
0206   TTree* m_tree = nullptr;
0207   ParticleType m_particleType;
0208 
0209   std::set<std::string> m_seenCandidateKeys;
0210 
0211   long long m_numberOfUniqueCandidates = 0;
0212   long long m_numberOfDuplicateCandidates = 0;
0213   bool m_currentEntryIsUnique = false;
0214 
0215   bool m_useDeltaBCOCut = false;
0216   long long m_minDeltaBCO = 0;
0217   long long m_maxDeltaBCO = 350;
0218   long long m_deltaBCO = -1;
0219   bool m_passesDeltaBCO = true;
0220 
0221   int m_runNumber = -1;
0222   int m_eventNumber = -1;
0223   Long64_t m_bco = -1;
0224   Long64_t m_eventBCO = -1;
0225 
0226   float m_mass = 0.0;
0227   float m_pt = 0.0;
0228   float m_eta = 0.0;
0229   float m_phi = 0.0;
0230   float m_decayLength = 0.0;
0231   float m_dira = 0.0;
0232   float m_ip = 0.0;
0233   float m_x = 0.0;
0234   float m_y = 0.0;
0235   float m_z = 0.0;
0236   float m_chi2 = 0.0;
0237   float m_rapidity = 0.0;
0238 };
0239 
0240 #endif