CMS 3D CMS Logo

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