CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HGCalTriggerNtupleHGCClusters.cc
Go to the documentation of this file.
1 
2 #include <algorithm>
8 
10 public:
13  void initialize(TTree&, const edm::ParameterSet&, edm::ConsumesCollector&&) final;
14  void fill(const edm::Event& e, const HGCalTriggerNtupleEventSetup& es) final;
15 
16 private:
17  void clear() final;
18 
22 
23  int cl_n_;
24  std::vector<uint32_t> cl_id_;
25  std::vector<float> cl_mipPt_;
26  std::vector<float> cl_pt_;
27  std::vector<float> cl_energy_;
28  std::vector<float> cl_eta_;
29  std::vector<float> cl_phi_;
30  std::vector<int> cl_layer_;
31  std::vector<int> cl_subdet_;
32  std::vector<int> cl_cells_n_;
33  std::vector<std::vector<uint32_t>> cl_cells_id_;
34  std::vector<uint32_t> cl_multicluster_id_;
36 };
37 
39 
40 HGCalTriggerNtupleHGCClusters::HGCalTriggerNtupleHGCClusters(const edm::ParameterSet& conf)
41  : HGCalTriggerNtupleBase(conf),
42  filter_clusters_in_multiclusters_(conf.getParameter<bool>("FilterClustersInMulticlusters")) {
43  accessEventSetup_ = false;
44 }
45 
47  const edm::ParameterSet& conf,
48  edm::ConsumesCollector&& collector) {
49  clusters_token_ = collector.consumes<l1t::HGCalClusterBxCollection>(conf.getParameter<edm::InputTag>("Clusters"));
51  collector.consumes<l1t::HGCalMulticlusterBxCollection>(conf.getParameter<edm::InputTag>("Multiclusters"));
52 
53  std::string prefix(conf.getUntrackedParameter<std::string>("Prefix", "cl"));
54 
55  std::string bname;
56  auto withPrefix([&prefix, &bname](char const* vname) -> char const* {
57  bname = prefix + "_" + vname;
58  return bname.c_str();
59  });
60 
61  // note: can't use withPrefix() twice within a same statement because bname gets overwritten
62  tree.Branch(withPrefix("n"), &cl_n_, (prefix + "_n/I").c_str());
63  tree.Branch(withPrefix("id"), &cl_id_);
64  tree.Branch(withPrefix("mipPt"), &cl_mipPt_);
65  tree.Branch(withPrefix("pt"), &cl_pt_);
66  tree.Branch(withPrefix("energy"), &cl_energy_);
67  tree.Branch(withPrefix("eta"), &cl_eta_);
68  tree.Branch(withPrefix("phi"), &cl_phi_);
69  tree.Branch(withPrefix("layer"), &cl_layer_);
70  tree.Branch(withPrefix("subdet"), &cl_subdet_);
71  tree.Branch(withPrefix("cells_n"), &cl_cells_n_);
72  tree.Branch(withPrefix("cells_id"), &cl_cells_id_);
73  tree.Branch(withPrefix("multicluster_id"), &cl_multicluster_id_);
74  tree.Branch(withPrefix("multicluster_pt"), &cl_multicluster_pt_);
75 }
76 
78  // retrieve clusters
80  e.getByToken(clusters_token_, clusters_h);
81  const l1t::HGCalClusterBxCollection& clusters = *clusters_h;
83  e.getByToken(multiclusters_token_, multiclusters_h);
84  const l1t::HGCalMulticlusterBxCollection& multiclusters = *multiclusters_h;
85 
87 
88  // Associate cells to clusters
89  std::unordered_map<uint32_t, l1t::HGCalMulticlusterBxCollection::const_iterator> cluster2multicluster;
90  for (auto mcl_itr = multiclusters.begin(0); mcl_itr != multiclusters.end(0); mcl_itr++) {
91  // loop on 2D clusters inside 3D clusters
92  for (const auto& cl_ptr : mcl_itr->constituents()) {
93  cluster2multicluster.emplace(cl_ptr.second->detId(), mcl_itr);
94  }
95  }
96 
97  clear();
98  for (auto cl_itr = clusters.begin(0); cl_itr != clusters.end(0); cl_itr++) {
99  auto mcl_itr = cluster2multicluster.find(cl_itr->detId());
100  uint32_t mcl_id = (mcl_itr != cluster2multicluster.end() ? mcl_itr->second->detId() : 0);
101  float mcl_pt = (mcl_itr != cluster2multicluster.end() ? mcl_itr->second->pt() : 0.);
102  if (filter_clusters_in_multiclusters_ && mcl_id == 0)
103  continue;
104  cl_n_++;
105  cl_mipPt_.emplace_back(cl_itr->mipPt());
106  // physical values
107  cl_pt_.emplace_back(cl_itr->pt());
108  cl_energy_.emplace_back(cl_itr->energy());
109  cl_eta_.emplace_back(cl_itr->eta());
110  cl_phi_.emplace_back(cl_itr->phi());
111 
112  cl_id_.emplace_back(cl_itr->detId());
113  cl_layer_.emplace_back(triggerTools_.layerWithOffset(cl_itr->detId()));
114  cl_subdet_.emplace_back(cl_itr->subdetId());
115  cl_cells_n_.emplace_back(cl_itr->constituents().size());
116  // Retrieve indices of trigger cells inside cluster
117  cl_cells_id_.emplace_back(cl_itr->constituents().size());
119  cl_itr->constituents_begin(),
120  cl_itr->constituents_end(),
121  cl_cells_id_.back().begin(),
122  [](const std::pair<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& id_tc) { return id_tc.second->detId(); });
123  cl_multicluster_id_.emplace_back(mcl_id);
124  cl_multicluster_pt_.emplace_back(mcl_pt);
125  }
126 }
127 
129  cl_n_ = 0;
130  cl_id_.clear();
131  cl_mipPt_.clear();
132  cl_pt_.clear();
133  cl_energy_.clear();
134  cl_eta_.clear();
135  cl_phi_.clear();
136  cl_layer_.clear();
137  cl_subdet_.clear();
138  cl_cells_n_.clear();
139  cl_cells_id_.clear();
140  cl_multicluster_id_.clear();
141  cl_multicluster_pt_.clear();
142 }
const_iterator end(int bx) const
T getUntrackedParameter(std::string const &, T const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
void fill(const edm::Event &e, const HGCalTriggerNtupleEventSetup &es) final
unsigned layerWithOffset(const DetId &) const
edm::ESHandle< HGCalTriggerGeometryBase > geometry
std::vector< std::vector< uint32_t > > cl_cells_id_
void setGeometry(const HGCalTriggerGeometryBase *const)
void initialize(TTree &, const edm::ParameterSet &, edm::ConsumesCollector &&) final
T const * product() const
Definition: ESHandle.h:86
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
HGCalTriggerNtupleHGCClusters(const edm::ParameterSet &conf)
#define DEFINE_EDM_PLUGIN(factory, type, name)
const_iterator begin(int bx) const
unsigned transform(const HcalDetId &id, unsigned transformCode)