Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2021-12-13 13:14:09

0001 #include "PhysicsTools/NanoAOD/plugins/LumiOutputBranches.h"
0002 
0003 #include <iostream>
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 LumiOutputBranches::defineBranchesFromFirstEvent(const nanoaod::FlatTable &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     switch (tab.columnType(i)) {
0016       case nanoaod::FlatTable::ColumnType::Float:
0017         m_floatBranches.emplace_back(var, tab.columnDoc(i), "F");
0018         break;
0019       case nanoaod::FlatTable::ColumnType::Int:
0020         m_intBranches.emplace_back(var, tab.columnDoc(i), "I");
0021         break;
0022       case nanoaod::FlatTable::ColumnType::UInt8:
0023         m_uint8Branches.emplace_back(var, tab.columnDoc(i), "b");
0024         break;
0025       case nanoaod::FlatTable::ColumnType::Bool:
0026         m_uint8Branches.emplace_back(var, tab.columnDoc(i), "O");
0027         break;
0028       default:
0029         throw cms::Exception("LogicError", "Unsupported type");
0030     }
0031   }
0032 }
0033 
0034 void LumiOutputBranches::branch(TTree &tree) {
0035   if (!m_singleton) {
0036     if (m_extension == IsExtension) {
0037       m_counterBranch = tree.FindBranch(("n" + m_baseName).c_str());
0038       if (!m_counterBranch) {
0039         throw cms::Exception("LogicError",
0040                              "Trying to save an extension table for " + m_baseName +
0041                                  " before having saved the corresponding main table\n");
0042       }
0043     } else {
0044       if (tree.FindBranch(("n" + m_baseName).c_str()) != nullptr) {
0045         throw cms::Exception("LogicError", "Trying to save multiple main tables for " + m_baseName + "\n");
0046       }
0047       m_counterBranch = tree.Branch(("n" + m_baseName).c_str(), &m_counter, ("n" + m_baseName + "/i").c_str());
0048       m_counterBranch->SetTitle(m_doc.c_str());
0049     }
0050   }
0051   std::string varsize = m_singleton ? "" : "[n" + m_baseName + "]";
0052   for (std::vector<NamedBranchPtr> *branches : {&m_floatBranches, &m_intBranches, &m_uint8Branches}) {
0053     for (auto &pair : *branches) {
0054       std::string branchName = makeBranchName(m_baseName, pair.name);
0055       pair.branch =
0056           tree.Branch(branchName.c_str(), (void *)nullptr, (branchName + varsize + "/" + pair.rootTypeCode).c_str());
0057       pair.branch->SetTitle(pair.title.c_str());
0058     }
0059   }
0060 }
0061 
0062 void LumiOutputBranches::fill(const edm::LuminosityBlockForOutput &iLumi, TTree &tree, bool extensions) {
0063   if (m_extension != DontKnowYetIfMainOrExtension) {
0064     if (extensions != m_extension)
0065       return;  // do nothing, wait to be called with the proper flag
0066   }
0067 
0068   edm::Handle<nanoaod::FlatTable> handle;
0069   iLumi.getByToken(m_token, handle);
0070   const nanoaod::FlatTable &tab = *handle;
0071   m_counter = tab.size();
0072   m_singleton = tab.singleton();
0073   if (!m_branchesBooked) {
0074     m_extension = tab.extension() ? IsExtension : IsMain;
0075     if (extensions != m_extension)
0076       return;  // do nothing, wait to be called with the proper flag
0077     defineBranchesFromFirstEvent(tab);
0078     m_doc = tab.doc();
0079     m_branchesBooked = true;
0080     branch(tree);
0081   }
0082   if (!m_singleton && m_extension == IsExtension) {
0083     if (m_counter != *reinterpret_cast<UInt_t *>(m_counterBranch->GetAddress())) {
0084       throw cms::Exception("LogicError",
0085                            "Mismatch in number of entries between extension and main table for " + tab.name());
0086     }
0087   }
0088   for (auto &pair : m_floatBranches)
0089     fillColumn<float>(pair, tab);
0090   for (auto &pair : m_intBranches)
0091     fillColumn<int>(pair, tab);
0092   for (auto &pair : m_uint8Branches)
0093     fillColumn<uint8_t>(pair, tab);
0094 }