Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-17 09:21:20

0001 /// ===========================================================================
0002 /*! \file   TrksInJetQABaseManager.cc
0003  *  \author Derek Anderson
0004  *  \date   04.03.2024
0005  *
0006  *  Base hist manager submodule for the TrksInJetQA module which
0007  *  consolidates methods/data common to all of the hist managers
0008  */
0009 /// ===========================================================================
0010 
0011 #include "TrksInJetQABaseManager.h"
0012 
0013 // ctor/dtor ==================================================================
0014 
0015 // ----------------------------------------------------------------------------
0016 //! Default ctor
0017 // ----------------------------------------------------------------------------
0018 TrksInJetQABaseManager::TrksInJetQABaseManager(TrksInJetQAConfig& config,
0019                                                TrksInJetQAHist& hist)
0020   : m_config(config)
0021   , m_hist(hist)
0022 {};
0023 
0024 // public methods =============================================================
0025 
0026 // ----------------------------------------------------------------------------
0027 //! Define and generate histograms for manager
0028 // ----------------------------------------------------------------------------
0029 void TrksInJetQABaseManager::MakeHistograms(const std::string& prefix,
0030                                             const std::string& suffix)
0031 {
0032   DefineHistograms();
0033   BuildHistograms(prefix, suffix);
0034 }  // end 'MakeHistograms(std::string)'
0035 
0036 // ----------------------------------------------------------------------------
0037 //! Save histograms to output file
0038 // ----------------------------------------------------------------------------
0039 /*! Note that this only relevant if output
0040  *  mode is OutMode::File.
0041  */
0042 void TrksInJetQABaseManager::SaveHistograms(TDirectory* topDir, const std::string& outDirName)
0043 {
0044   TDirectory* outDir = topDir->mkdir(outDirName.c_str());
0045   if (!outDir)
0046   {
0047     std::cerr << PHWHERE << ": PANIC: unable to make output directory!" << std::endl;
0048     exit(1);
0049   }
0050 
0051   outDir->cd();
0052   for (const auto& hist1D : m_mapHist1D)
0053   {
0054     hist1D.second->Write();
0055   }
0056   for (const auto& hist2D : m_mapHist2D)
0057   {
0058     hist2D.second->Write();
0059   }
0060 }  // end 'SaveHistograms(TDirectory*, std::string)'
0061 
0062 // ----------------------------------------------------------------------------
0063 //! Grab histograms from manager
0064 // ----------------------------------------------------------------------------
0065 void TrksInJetQABaseManager::GrabHistograms(std::vector<TH1D*>& vecOutHist1D,
0066                                             std::vector<TH2D*>& vecOutHist2D)
0067 {
0068   for (const auto& hist1D : m_mapHist1D)
0069   {
0070     vecOutHist1D.push_back(hist1D.second);
0071   }
0072   for (const auto& hist2D : m_mapHist2D)
0073   {
0074     vecOutHist2D.push_back(hist2D.second);
0075   }
0076 }  // end 'GrabHistograms(std::vector<TH1D*>&, std::vector<TH2D*>&)'
0077 
0078 // private methods ============================================================
0079 
0080 // ----------------------------------------------------------------------------
0081 //! Build histograms from definitions
0082 // ----------------------------------------------------------------------------
0083 /*! Note that the specific histogram definitions are
0084  *  implemented in the derived classes.
0085  */
0086 void TrksInJetQABaseManager::BuildHistograms(const std::string& prefix,
0087                                              const std::string& suffix)
0088 {
0089   // build 1d histograms
0090   for (const auto& histType : m_mapHistTypes)
0091   {
0092     for (const auto& histDef1D : m_mapHistDef1D)
0093     {
0094       // make name
0095       std::string sHistName(prefix + "_");
0096       sHistName += histType.second;
0097       sHistName += std::get<0>(histDef1D.second);
0098       sHistName += "_";
0099       sHistName += suffix;
0100 
0101       // make sure histogram name is lower case
0102       std::transform(sHistName.begin(),
0103                      sHistName.end(),
0104                      sHistName.begin(),
0105                      ::tolower);
0106       std::regex_replace(
0107           sHistName,
0108           std::regex("__"),
0109           "_");
0110 
0111       // create histogram
0112       m_mapHist1D.emplace(Index(histType.first, histDef1D.first),
0113                           new TH1D(sHistName.data(),
0114                                    "",
0115                                    std::get<1>(histDef1D.second).first,
0116                                    std::get<1>(histDef1D.second).second.first,
0117                                    std::get<1>(histDef1D.second).second.second));
0118     }  // end 1D hist loop
0119 
0120     // build 2d histograms
0121     for (const auto& histDef2D : m_mapHistDef2D)
0122     {
0123       // make name
0124       std::string sHistName(prefix + "_");
0125       sHistName += histType.second;
0126       sHistName += std::get<0>(histDef2D.second);
0127       sHistName += "_";
0128       sHistName += suffix;
0129 
0130       // make sure histogram name is lower case
0131       std::transform(sHistName.begin(),
0132                      sHistName.end(),
0133                      sHistName.begin(),
0134                      ::tolower);
0135       std::regex_replace(
0136           sHistName,
0137           std::regex("__"),
0138           "_");
0139 
0140       // create histogram
0141       m_mapHist2D.emplace(Index(histType.first, histDef2D.first),
0142                           new TH2D(sHistName.data(),
0143                                    "",
0144                                    std::get<1>(histDef2D.second).first,
0145                                    std::get<1>(histDef2D.second).second.first,
0146                                    std::get<1>(histDef2D.second).second.second,
0147                                    std::get<2>(histDef2D.second).first,
0148                                    std::get<2>(histDef2D.second).second.first,
0149                                    std::get<2>(histDef2D.second).second.second));
0150     }  // end hist loop
0151   }  // end type loop
0152 }  // end 'BuildHistograms(std::string)'
0153 
0154 // private helper methods =====================================================
0155 
0156 // ----------------------------------------------------------------------------
0157 //! Check if a layer is in the MVTX
0158 // ----------------------------------------------------------------------------
0159 bool TrksInJetQABaseManager::IsInMvtx(const uint16_t layer) const
0160 {
0161   return (layer < m_config.nMvtxLayer);
0162 }  // end 'IsInMvtx(uint16_t)'
0163 
0164 // ----------------------------------------------------------------------------
0165 //! Check if a layer is in the INTT
0166 // ----------------------------------------------------------------------------
0167 bool TrksInJetQABaseManager::IsInIntt(const uint16_t layer) const
0168 {
0169   return ((layer >= m_config.nMvtxLayer) &&
0170           (layer < m_config.nInttLayer + m_config.nMvtxLayer));
0171 }  // end 'IsInIntt(uint16_t)'
0172 
0173 // ----------------------------------------------------------------------------
0174 //! Check if a layer is in the TPC
0175 // ----------------------------------------------------------------------------
0176 bool TrksInJetQABaseManager::IsInTpc(const uint16_t layer) const
0177 {
0178   return (layer >= m_config.nMvtxLayer + m_config.nInttLayer);
0179 }  // end 'IsInTpc(uint16_t)'
0180 
0181 // ----------------------------------------------------------------------------
0182 //! Create an index based on object type and histogram index
0183 // ----------------------------------------------------------------------------
0184 int TrksInJetQABaseManager::Index(const int type, const int hist) const
0185 {
0186   return (100 * type) + hist;
0187 }  // end 'Index(int, int)'
0188 
0189 // end ========================================================================