CMS 3D CMS Logo

HGCalTriggerNtupleHGCClusters.cc
Go to the documentation of this file.
1 
2 #include <algorithm>
9 
11 public:
14  void initialize(TTree&, const edm::ParameterSet&, edm::ConsumesCollector&&) final;
15  void fill(const edm::Event& e, const edm::EventSetup& es) final;
16 
17 private:
18  void clear() final;
19 
23 
24  int cl_n_;
25  std::vector<uint32_t> cl_id_;
26  std::vector<float> cl_mipPt_;
27  std::vector<float> cl_pt_;
28  std::vector<float> cl_energy_;
29  std::vector<float> cl_eta_;
30  std::vector<float> cl_phi_;
31  std::vector<int> cl_layer_;
32  std::vector<int> cl_subdet_;
33  std::vector<int> cl_cells_n_;
34  std::vector<std::vector<uint32_t>> cl_cells_id_;
35  std::vector<uint32_t> cl_multicluster_id_;
36  std::vector<float> cl_multicluster_pt_;
37 };
38 
40 
42  : HGCalTriggerNtupleBase(conf),
43  filter_clusters_in_multiclusters_(conf.getParameter<bool>("FilterClustersInMulticlusters")) {}
44 
46  const edm::ParameterSet& conf,
47  edm::ConsumesCollector&& collector) {
48  clusters_token_ = collector.consumes<l1t::HGCalClusterBxCollection>(conf.getParameter<edm::InputTag>("Clusters"));
50  collector.consumes<l1t::HGCalMulticlusterBxCollection>(conf.getParameter<edm::InputTag>("Multiclusters"));
51 
52  tree.Branch("cl_n", &cl_n_, "cl_n/I");
53  tree.Branch("cl_id", &cl_id_);
54  tree.Branch("cl_mipPt", &cl_mipPt_);
55  tree.Branch("cl_pt", &cl_pt_);
56  tree.Branch("cl_energy", &cl_energy_);
57  tree.Branch("cl_eta", &cl_eta_);
58  tree.Branch("cl_phi", &cl_phi_);
59  tree.Branch("cl_layer", &cl_layer_);
60  tree.Branch("cl_subdet", &cl_subdet_);
61  tree.Branch("cl_cells_n", &cl_cells_n_);
62  tree.Branch("cl_cells_id", &cl_cells_id_);
63  tree.Branch("cl_multicluster_id", &cl_multicluster_id_);
64  tree.Branch("cl_multicluster_pt", &cl_multicluster_pt_);
65 }
66 
68  // retrieve clusters
70  e.getByToken(clusters_token_, clusters_h);
71  const l1t::HGCalClusterBxCollection& clusters = *clusters_h;
73  e.getByToken(multiclusters_token_, multiclusters_h);
74  const l1t::HGCalMulticlusterBxCollection& multiclusters = *multiclusters_h;
75 
76  // retrieve geometry
78  es.get<CaloGeometryRecord>().get(geometry);
79 
81 
82  // Associate cells to clusters
83  std::unordered_map<uint32_t, l1t::HGCalMulticlusterBxCollection::const_iterator> cluster2multicluster;
84  for (auto mcl_itr = multiclusters.begin(0); mcl_itr != multiclusters.end(0); mcl_itr++) {
85  // loop on 2D clusters inside 3D clusters
86  for (const auto& cl_ptr : mcl_itr->constituents()) {
87  cluster2multicluster.emplace(cl_ptr.second->detId(), mcl_itr);
88  }
89  }
90 
91  clear();
92  for (auto cl_itr = clusters.begin(0); cl_itr != clusters.end(0); cl_itr++) {
93  auto mcl_itr = cluster2multicluster.find(cl_itr->detId());
94  uint32_t mcl_id = (mcl_itr != cluster2multicluster.end() ? mcl_itr->second->detId() : 0);
95  float mcl_pt = (mcl_itr != cluster2multicluster.end() ? mcl_itr->second->pt() : 0.);
96  if (filter_clusters_in_multiclusters_ && mcl_id == 0)
97  continue;
98  cl_n_++;
99  cl_mipPt_.emplace_back(cl_itr->mipPt());
100  // physical values
101  cl_pt_.emplace_back(cl_itr->pt());
102  cl_energy_.emplace_back(cl_itr->energy());
103  cl_eta_.emplace_back(cl_itr->eta());
104  cl_phi_.emplace_back(cl_itr->phi());
105 
106  cl_id_.emplace_back(cl_itr->detId());
107  cl_layer_.emplace_back(triggerTools_.layerWithOffset(cl_itr->detId()));
108  cl_subdet_.emplace_back(cl_itr->subdetId());
109  cl_cells_n_.emplace_back(cl_itr->constituents().size());
110  // Retrieve indices of trigger cells inside cluster
111  cl_cells_id_.emplace_back(cl_itr->constituents().size());
113  cl_itr->constituents_begin(),
114  cl_itr->constituents_end(),
115  cl_cells_id_.back().begin(),
116  [](const std::pair<uint32_t, edm::Ptr<l1t::HGCalTriggerCell>>& id_tc) { return id_tc.second->detId(); });
117  cl_multicluster_id_.emplace_back(mcl_id);
118  cl_multicluster_pt_.emplace_back(mcl_pt);
119  }
120 }
121 
123  cl_n_ = 0;
124  cl_id_.clear();
125  cl_mipPt_.clear();
126  cl_pt_.clear();
127  cl_energy_.clear();
128  cl_eta_.clear();
129  cl_phi_.clear();
130  cl_layer_.clear();
131  cl_subdet_.clear();
132  cl_cells_n_.clear();
133  cl_cells_id_.clear();
134  cl_multicluster_id_.clear();
135  cl_multicluster_pt_.clear();
136 }
const_iterator end(int bx) const
T getParameter(std::string const &) const
void eventSetup(const edm::EventSetup &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
unsigned layerWithOffset(const DetId &) const
std::vector< std::vector< uint32_t > > cl_cells_id_
void initialize(TTree &, const edm::ParameterSet &, edm::ConsumesCollector &&) final
ESHandle< TrackerGeometry > geometry
T get() const
Definition: EventSetup.h:71
HGCalTriggerNtupleHGCClusters(const edm::ParameterSet &conf)
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: tree.py:1
const_iterator begin(int bx) const
void fill(const edm::Event &e, const edm::EventSetup &es) final