CMS 3D CMS Logo

HGCalTriggerNtupleHGCMulticlusters.cc
Go to the documentation of this file.
6 
8 public:
11  void initialize(TTree&, const edm::ParameterSet&, edm::ConsumesCollector&&) final;
12  void fill(const edm::Event& e, const edm::EventSetup& es) final;
13 
14 private:
15  void clear() final;
16 
18 
19  std::unique_ptr<HGCalTriggerClusterIdentificationBase> id_;
20 
21  int cl3d_n_;
22  std::vector<uint32_t> cl3d_id_;
23  std::vector<float> cl3d_pt_;
24  std::vector<float> cl3d_energy_;
25  std::vector<float> cl3d_eta_;
26  std::vector<float> cl3d_phi_;
27  std::vector<int> cl3d_clusters_n_;
28  std::vector<std::vector<uint32_t>> cl3d_clusters_id_;
29  // cluster shower shapes
30  std::vector<int> cl3d_showerlength_;
31  std::vector<int> cl3d_coreshowerlength_;
32  std::vector<int> cl3d_firstlayer_;
33  std::vector<int> cl3d_maxlayer_;
34  std::vector<float> cl3d_seetot_;
35  std::vector<float> cl3d_seemax_;
36  std::vector<float> cl3d_spptot_;
37  std::vector<float> cl3d_sppmax_;
38  std::vector<float> cl3d_szz_;
39  std::vector<float> cl3d_srrtot_;
40  std::vector<float> cl3d_srrmax_;
41  std::vector<float> cl3d_srrmean_;
42  std::vector<float> cl3d_emaxe_;
43  std::vector<float> cl3d_bdteg_;
44  std::vector<int> cl3d_quality_;
45 };
46 
48 
50  : HGCalTriggerNtupleBase(conf) {}
51 
53  const edm::ParameterSet& conf,
54  edm::ConsumesCollector&& collector) {
56  collector.consumes<l1t::HGCalMulticlusterBxCollection>(conf.getParameter<edm::InputTag>("Multiclusters"));
57  id_ = std::unique_ptr<HGCalTriggerClusterIdentificationBase>{
58  HGCalTriggerClusterIdentificationFactory::get()->create("HGCalTriggerClusterIdentificationBDT")};
59  id_->initialize(conf.getParameter<edm::ParameterSet>("EGIdentification"));
60 
61  tree.Branch("cl3d_n", &cl3d_n_, "cl3d_n/I");
62  tree.Branch("cl3d_id", &cl3d_id_);
63  tree.Branch("cl3d_pt", &cl3d_pt_);
64  tree.Branch("cl3d_energy", &cl3d_energy_);
65  tree.Branch("cl3d_eta", &cl3d_eta_);
66  tree.Branch("cl3d_phi", &cl3d_phi_);
67  tree.Branch("cl3d_clusters_n", &cl3d_clusters_n_);
68  tree.Branch("cl3d_clusters_id", &cl3d_clusters_id_);
69  tree.Branch("cl3d_showerlength", &cl3d_showerlength_);
70  tree.Branch("cl3d_coreshowerlength", &cl3d_coreshowerlength_);
71  tree.Branch("cl3d_firstlayer", &cl3d_firstlayer_);
72  tree.Branch("cl3d_maxlayer", &cl3d_maxlayer_);
73  tree.Branch("cl3d_seetot", &cl3d_seetot_);
74  tree.Branch("cl3d_seemax", &cl3d_seemax_);
75  tree.Branch("cl3d_spptot", &cl3d_spptot_);
76  tree.Branch("cl3d_sppmax", &cl3d_sppmax_);
77  tree.Branch("cl3d_szz", &cl3d_szz_);
78  tree.Branch("cl3d_srrtot", &cl3d_srrtot_);
79  tree.Branch("cl3d_srrmax", &cl3d_srrmax_);
80  tree.Branch("cl3d_srrmean", &cl3d_srrmean_);
81  tree.Branch("cl3d_emaxe", &cl3d_emaxe_);
82  tree.Branch("cl3d_bdteg", &cl3d_bdteg_);
83  tree.Branch("cl3d_quality", &cl3d_quality_);
84 }
85 
87  // retrieve clusters 3D
89  e.getByToken(multiclusters_token_, multiclusters_h);
90  const l1t::HGCalMulticlusterBxCollection& multiclusters = *multiclusters_h;
91 
92  // retrieve geometry
94  es.get<CaloGeometryRecord>().get(geometry);
95 
96  clear();
97  for (auto cl3d_itr = multiclusters.begin(0); cl3d_itr != multiclusters.end(0); cl3d_itr++) {
98  cl3d_n_++;
99  cl3d_id_.emplace_back(cl3d_itr->detId());
100  // physical values
101  cl3d_pt_.emplace_back(cl3d_itr->pt());
102  cl3d_energy_.emplace_back(cl3d_itr->energy());
103  cl3d_eta_.emplace_back(cl3d_itr->eta());
104  cl3d_phi_.emplace_back(cl3d_itr->phi());
105  cl3d_clusters_n_.emplace_back(cl3d_itr->constituents().size());
106  cl3d_showerlength_.emplace_back(cl3d_itr->showerLength());
107  cl3d_coreshowerlength_.emplace_back(cl3d_itr->coreShowerLength());
108  cl3d_firstlayer_.emplace_back(cl3d_itr->firstLayer());
109  cl3d_maxlayer_.emplace_back(cl3d_itr->maxLayer());
110  cl3d_seetot_.emplace_back(cl3d_itr->sigmaEtaEtaTot());
111  cl3d_seemax_.emplace_back(cl3d_itr->sigmaEtaEtaMax());
112  cl3d_spptot_.emplace_back(cl3d_itr->sigmaPhiPhiTot());
113  cl3d_sppmax_.emplace_back(cl3d_itr->sigmaPhiPhiMax());
114  cl3d_szz_.emplace_back(cl3d_itr->sigmaZZ());
115  cl3d_srrtot_.emplace_back(cl3d_itr->sigmaRRTot());
116  cl3d_srrmax_.emplace_back(cl3d_itr->sigmaRRMax());
117  cl3d_srrmean_.emplace_back(cl3d_itr->sigmaRRMean());
118  cl3d_emaxe_.emplace_back(cl3d_itr->eMax() / cl3d_itr->energy());
119  cl3d_bdteg_.emplace_back(id_->value(*cl3d_itr));
120  cl3d_quality_.emplace_back(cl3d_itr->hwQual());
121 
122  // Retrieve indices of trigger cells inside cluster
123  cl3d_clusters_id_.emplace_back(cl3d_itr->constituents().size());
124  std::transform(cl3d_itr->constituents_begin(),
125  cl3d_itr->constituents_end(),
126  cl3d_clusters_id_.back().begin(),
127  [](const std::pair<uint32_t, edm::Ptr<l1t::HGCalCluster>>& id_cl) { return id_cl.second->detId(); });
128  }
129 }
130 
132  cl3d_n_ = 0;
133  cl3d_id_.clear();
134  cl3d_pt_.clear();
135  cl3d_energy_.clear();
136  cl3d_eta_.clear();
137  cl3d_phi_.clear();
138  cl3d_clusters_n_.clear();
139  cl3d_clusters_id_.clear();
140  cl3d_showerlength_.clear();
141  cl3d_coreshowerlength_.clear();
142  cl3d_firstlayer_.clear();
143  cl3d_maxlayer_.clear();
144  cl3d_seetot_.clear();
145  cl3d_seemax_.clear();
146  cl3d_spptot_.clear();
147  cl3d_sppmax_.clear();
148  cl3d_szz_.clear();
149  cl3d_srrtot_.clear();
150  cl3d_srrmax_.clear();
151  cl3d_srrmean_.clear();
152  cl3d_emaxe_.clear();
153  cl3d_bdteg_.clear();
154  cl3d_quality_.clear();
155 }
const_iterator end(int bx) const
T getParameter(std::string const &) const
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
HGCalTriggerNtupleHGCMulticlusters(const edm::ParameterSet &conf)
std::vector< std::vector< uint32_t > > cl3d_clusters_id_
std::unique_ptr< HGCalTriggerClusterIdentificationBase > id_
ESHandle< TrackerGeometry > geometry
void initialize(TTree &, const edm::ParameterSet &, edm::ConsumesCollector &&) final
void fill(const edm::Event &e, const edm::EventSetup &es) final
T get() const
Definition: EventSetup.h:71
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: tree.py:1
const_iterator begin(int bx) const
T get(const Candidate &c)
Definition: component.h:55