Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-08-12 02:29:05

0001 #include "L1TriggerScouting/Utilities/plugins/OrbitTableOutputBranches.h"
0002 
0003 #include <limits>
0004 
0005 namespace {
0006   std::string makeBranchName(const std::string &baseName, const std::string &leafName) {
0007     return baseName.empty() ? leafName : (leafName.empty() ? baseName : baseName + "_" + leafName);
0008   }
0009 }  // namespace
0010 
0011 void OrbitTableOutputBranches::defineBranchesFromFirstEvent(const l1ScoutingRun3::OrbitFlatTable &tab) {
0012   m_baseName = tab.name();
0013   for (size_t i = 0; i < tab.nColumns(); i++) {
0014     const std::string &var = tab.columnName(i);
0015     // get columnIndex
0016     int columnIndex = tab.columnIndex(var);
0017     if (columnIndex == -1)
0018       throw cms::Exception("LogicError", "Missing column in input for " + m_baseName + "_" + var);
0019 
0020     switch (tab.columnType(i)) {
0021       case l1ScoutingRun3::OrbitFlatTable::ColumnType::UInt8:
0022         m_uint8Branches.emplace_back(var, tab.columnDoc(i), "b", columnIndex);
0023         break;
0024       case l1ScoutingRun3::OrbitFlatTable::ColumnType::Int16:
0025         m_int16Branches.emplace_back(var, tab.columnDoc(i), "S", columnIndex);
0026         break;
0027       case l1ScoutingRun3::OrbitFlatTable::ColumnType::UInt16:
0028         m_uint16Branches.emplace_back(var, tab.columnDoc(i), "s", columnIndex);
0029         break;
0030       case l1ScoutingRun3::OrbitFlatTable::ColumnType::Int32:
0031         m_int32Branches.emplace_back(var, tab.columnDoc(i), "I", columnIndex);
0032         break;
0033       case l1ScoutingRun3::OrbitFlatTable::ColumnType::UInt32:
0034         m_uint32Branches.emplace_back(var, tab.columnDoc(i), "i", columnIndex);
0035         break;
0036       case l1ScoutingRun3::OrbitFlatTable::ColumnType::Bool:
0037         m_uint8Branches.emplace_back(var, tab.columnDoc(i), "O", columnIndex);
0038         break;
0039       case l1ScoutingRun3::OrbitFlatTable::ColumnType::Float:
0040         m_floatBranches.emplace_back(var, tab.columnDoc(i), "F", columnIndex);
0041         break;
0042       case l1ScoutingRun3::OrbitFlatTable::ColumnType::Double:
0043         m_doubleBranches.emplace_back(var, tab.columnDoc(i), "D", columnIndex);
0044         break;
0045       default:
0046         throw cms::Exception("LogicError", "Unsupported type");
0047     }
0048   }
0049 }
0050 
0051 void OrbitTableOutputBranches::branch(TTree &tree) {
0052   if (!m_singleton) {
0053     if (m_extension == IsExtension) {
0054       m_counterBranch = tree.FindBranch(("n" + m_baseName).c_str());
0055       if (!m_counterBranch) {
0056         throw cms::Exception("LogicError",
0057                              "Trying to save an extension table for " + m_baseName +
0058                                  " before having saved the corresponding main table\n");
0059       }
0060     } else {
0061       if (tree.FindBranch(("n" + m_baseName).c_str()) != nullptr) {
0062         throw cms::Exception("LogicError", "Trying to save multiple main tables for " + m_baseName + "\n");
0063       }
0064       m_counterBranch = tree.Branch(("n" + m_baseName).c_str(), &m_counter, ("n" + m_baseName + "/I").c_str());
0065       m_counterBranch->SetTitle(m_doc.c_str());
0066     }
0067   }
0068   std::string varsize = m_singleton ? "" : "[n" + m_baseName + "]";
0069   for (std::vector<NamedBranchPtr> *branches : {&m_uint8Branches,
0070                                                 &m_int16Branches,
0071                                                 &m_uint16Branches,
0072                                                 &m_int32Branches,
0073                                                 &m_uint32Branches,
0074                                                 &m_floatBranches,
0075                                                 &m_doubleBranches}) {
0076     for (auto &pair : *branches) {
0077       std::string branchName = makeBranchName(m_baseName, pair.name);
0078       pair.branch =
0079           tree.Branch(branchName.c_str(), (void *)nullptr, (branchName + varsize + "/" + pair.rootTypeCode).c_str());
0080       pair.branch->SetTitle(pair.title.c_str());
0081     }
0082   }
0083 }
0084 
0085 void OrbitTableOutputBranches::beginFill(const edm::OccurrenceForOutput &iWhatever, TTree &tree, bool extensions) {
0086   if (m_extension != DontKnowYetIfMainOrExtension) {
0087     if (extensions != m_extension)
0088       return;  // do nothing, wait to be called with the proper flag
0089   }
0090 
0091   iWhatever.getByToken(m_token, m_handle);
0092   m_table = m_handle.product();
0093   m_singleton = m_table->singleton();
0094   if (!m_branchesBooked) {
0095     m_extension = m_table->extension() ? IsExtension : IsMain;
0096     if (extensions != m_extension)
0097       return;  // do nothing, wait to be called with the proper flag
0098     defineBranchesFromFirstEvent(*m_table);
0099     m_doc = m_table->doc();
0100     m_branchesBooked = true;
0101     branch(tree);
0102   }
0103 }
0104 
0105 bool OrbitTableOutputBranches::hasBx(uint32_t bx) { return (m_table->size(bx) != 0); }
0106 
0107 void OrbitTableOutputBranches::fillBx(uint32_t bx) {
0108   m_counter = m_table->size(bx);
0109   if (m_counter != 0) {
0110     for (auto &pair : m_uint8Branches)
0111       fillColumn<uint8_t>(pair, bx);
0112     for (auto &pair : m_int16Branches)
0113       fillColumn<int16_t>(pair, bx);
0114     for (auto &pair : m_uint16Branches)
0115       fillColumn<uint16_t>(pair, bx);
0116     for (auto &pair : m_int32Branches)
0117       fillColumn<int32_t>(pair, bx);
0118     for (auto &pair : m_uint32Branches)
0119       fillColumn<uint32_t>(pair, bx);
0120     for (auto &pair : m_floatBranches)
0121       fillColumn<float>(pair, bx);
0122     for (auto &pair : m_doubleBranches)
0123       fillColumn<double>(pair, bx);
0124   }
0125 }
0126 
0127 void OrbitTableOutputBranches::endFill() { m_table = nullptr; }