CMS 3D CMS Logo

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