CMS 3D CMS Logo

ElectronStudy.cc
Go to the documentation of this file.
3 
16 
24 
29 
30 #include <TH1F.h>
31 
32 #include <algorithm>
33 #include <fstream>
34 #include <iostream>
35 #include <memory>
36 #include <string>
37 #include <vector>
38 
39 class ElectronStudy: public edm::one::EDAnalyzer<edm::one::SharedResources> {
40 
41 public:
42 
43  explicit ElectronStudy(const edm::ParameterSet& ps);
44  ~ElectronStudy() override {}
45 
46  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
47 
48 private:
49 
50  void analyze(edm::Event const&, edm::EventSetup const&) override;
51  void beginJob() override {}
52 
53  static const int NEtaBins_ = 3;
54  static const int NPBins_ = 8;
55  double pBins_[NPBins_+1], etaBins_[NEtaBins_+1];
56 
61 
63  bool histos_;
64  TH1F *histoR1[NPBins_+1][NEtaBins_+1];
65  TH1F *histoR2[NPBins_+1][NEtaBins_+1];
66  TH1F *histoR3[NPBins_+1][NEtaBins_+1];
67  TH1F *histoE1x1[NPBins_+1][NEtaBins_+1];
68  TH1F *histoE3x3[NPBins_+1][NEtaBins_+1];
69  TH1F *histoE5x5[NPBins_+1][NEtaBins_+1];
70 };
71 
73 
74  usesResource("TFileService");
75 
76  std::string g4Label = ps.getUntrackedParameter<std::string>("ModuleLabel","g4SimHits");
77  std::string hitLabEB= ps.getUntrackedParameter<std::string>("EBCollection","EcalHitsEB");
78  std::string hitLabEE= ps.getUntrackedParameter<std::string>("EECollection","EcalHitsEE");
79 
80 
81  tok_EBhit_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label,hitLabEB));
82  tok_EEhit_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label,hitLabEE));
83  tok_simTk_ = consumes<edm::SimTrackContainer>(edm::InputTag(g4Label));
84  tok_simVtx_ = consumes<edm::SimVertexContainer>(edm::InputTag(g4Label));
85 
86  hotZone_ = ps.getUntrackedParameter<int>("HotZone",0);
87  verbose_ = ps.getUntrackedParameter<int>("Verbosity",0);
88  edm::LogInfo("ElectronStudy") << "Module Label: " << g4Label << " Hits: "
89  << hitLabEB << ", " << hitLabEE;
90 
91  double tempP[NPBins_+1] = { 0.0, 10.0, 20.0, 40.0, 60.0,
92  100.0, 500.0, 1000.0, 10000.0};
93  double tempEta[NEtaBins_+1] = {0.0, 1.2, 1.6, 3.0};
94 
95  for (int i=0; i<NPBins_+1; i++) pBins_[i] = tempP[i];
96  for(int i=0; i<NEtaBins_+1; i++) etaBins_[i] = tempEta[i];
97 
99  if ( !tfile.isAvailable() ) {
100  edm::LogInfo("ElectronStudy") << "TFileService unavailable: no histograms";
101  histos_ = false;
102  } else {
103  char name[20], title[200], cpbin[30], cebin[30];
104  histos_ = true;
105  for (unsigned int i=0; i<NPBins_+1; ++i) {
106  if (i == 0) sprintf (cpbin, " All p");
107  else sprintf (cpbin, " p (%6.0f:%6.0f)", pBins_[i-1], pBins_[i]);
108  for (unsigned int j=0; j<NEtaBins_+1; ++j) {
109  if (j == 0) sprintf (cebin, " All #eta");
110  else sprintf (cebin, " #eta (%4.1f:%4.1f)", etaBins_[j-1], etaBins_[j]);
111  sprintf (name, "R1%d%d", i, j);
112  sprintf (title,"E1/E9 for %s%s", cpbin, cebin);
113  histoR1[i][j] = tfile->make<TH1F>(name, title, 100, 0., 2.);
114  histoR1[i][j]->GetXaxis()->SetTitle(title);
115  histoR1[i][j]->GetYaxis()->SetTitle("Tracks");
116  sprintf (name, "R2%d%d", i, j);
117  sprintf (title,"E1/E25 for %s%s", cpbin, cebin);
118  histoR2[i][j] = tfile->make<TH1F>(name, title, 100, 0., 2.);
119  histoR2[i][j]->GetXaxis()->SetTitle(title);
120  histoR2[i][j]->GetYaxis()->SetTitle("Tracks");
121  sprintf (name, "R3%d%d", i, j);
122  sprintf (title,"E9/E25 for %s%s", cpbin, cebin);
123  histoR3[i][j] = tfile->make<TH1F>(name, title, 100, 0., 2.);
124  histoR3[i][j]->GetXaxis()->SetTitle(title);
125  histoR3[i][j]->GetYaxis()->SetTitle("Tracks");
126  sprintf (name, "E1x1%d%d", i, j);
127  sprintf (title,"E1/P for %s%s", cpbin, cebin);
128  histoE1x1[i][j] = tfile->make<TH1F>(name, title, 100, 0., 2.);
129  histoE1x1[i][j]->GetXaxis()->SetTitle(title);
130  histoE1x1[i][j]->GetYaxis()->SetTitle("Tracks");
131  sprintf (name, "E3x3%d%d", i, j);
132  sprintf (title,"E9/P for %s%s", cpbin, cebin);
133  histoE3x3[i][j] = tfile->make<TH1F>(name, title, 100, 0., 2.);
134  histoE3x3[i][j]->GetXaxis()->SetTitle(title);
135  histoE3x3[i][j]->GetYaxis()->SetTitle("Tracks");
136  sprintf (name, "E5x5%d%d", i, j);
137  sprintf (title,"E25/P for %s%s", cpbin, cebin);
138  histoE5x5[i][j] = tfile->make<TH1F>(name, title, 100, 0., 2.);
139  histoE5x5[i][j]->GetXaxis()->SetTitle(title);
140  histoE5x5[i][j]->GetYaxis()->SetTitle("Tracks");
141  }
142  }
143  }
144 }
145 
147 
149  desc.addUntracked<std::string>("ModuleLabel","g4SimHits");
150  desc.addUntracked<std::string>("EBCollection","EcalHitsEB");
151  desc.addUntracked<std::string>("EECollection","EcalHitsEE");
152  desc.addUntracked<int>("Verbosity",0);
153  descriptions.add("electronStudy",desc);
154 }
155 
157 
158  if (verbose_ > 1)
159  edm::LogVerbatim("IsoTrack") << "Run = " << iEvent.id().run()
160  << " Event = " << iEvent.id().event();
161 
162  // get Geometry, B-field, Topology
164  iSetup.get<IdealMagneticFieldRecord>().get(bFieldH);
165  const MagneticField* bField = bFieldH.product();
166 
168  iSetup.get<CaloGeometryRecord>().get(pG);
169  const CaloGeometry* geo = pG.product();
170 
171  edm::ESHandle<CaloTopology> theCaloTopology;
172  iSetup.get<CaloTopologyRecord>().get(theCaloTopology);
173  const CaloTopology* caloTopology = theCaloTopology.product();
174 
175  // get PCaloHits for ecal barrel
177  iEvent.getByToken(tok_EBhit_,caloHitEB);
178 
179  // get PCaloHits for ecal endcap
181  iEvent.getByToken(tok_EEhit_,caloHitEE);
182 
183  // get sim tracks
185  iEvent.getByToken(tok_simTk_, SimTk);
186 
187  // get sim vertices
189  iEvent.getByToken(tok_simVtx_, SimVtx);
190 
191  if (verbose_>0)
192  edm::LogVerbatim("IsoTrack") << "ElectronStudy: hits valid[EB]: "
193  << caloHitEB.isValid() << " valid[EE]: "
194  << caloHitEE.isValid();
195 
196  if (caloHitEB.isValid() && caloHitEE.isValid()) {
197  unsigned int indx;
198  if (verbose_>2) {
199  edm::PCaloHitContainer::const_iterator ihit;
200  for (ihit=caloHitEB->begin(),indx=0; ihit!=caloHitEB->end();
201  ihit++, indx++) {
202  EBDetId id = ihit->id();
203  edm::LogVerbatim("IsoTrack") << "Hit[" << indx << "] " << id << " E "
204  << ihit->energy() << " T " << ihit->time();
205  }
206  for (ihit=caloHitEE->begin(),indx=0; ihit!=caloHitEE->end();
207  ihit++,indx++) {
208  EEDetId id = ihit->id();
209  edm::LogVerbatim("IsoTrack") << "Hit[" << indx << "] " << id << " E "
210  << ihit->energy() << " T " << ihit->time();
211  }
212  }
213  edm::SimTrackContainer::const_iterator simTrkItr=SimTk->begin();
214  for (indx=0; simTrkItr!= SimTk->end(); simTrkItr++,indx++) {
215  if (verbose_>0)
216  edm::LogVerbatim("IsoTrack") << "ElectronStudy: Track[" << indx
217  << "] ID " << simTrkItr->trackId()
218  << " type " << simTrkItr->type()
219  << " charge " << simTrkItr->charge()
220  << " p " << simTrkItr->momentum()
221  << " Generator Index "
222  << simTrkItr->genpartIndex() << " vertex "
223  << simTrkItr->vertIndex();
224  if (std::abs(simTrkItr->type()) == 11 && simTrkItr->vertIndex() != -1) {
225  int thisTrk = simTrkItr->trackId();
226  spr::propagatedTrackDirection trkD = spr::propagateCALO(thisTrk, SimTk, SimVtx, geo, bField, (verbose_>1));
227  if (trkD.okECAL) {
228  const DetId isoCell = trkD.detIdECAL;
229  DetId hotCell = isoCell;
230  if (hotZone_ > 0) hotCell = spr::hotCrystal(isoCell, caloHitEB, caloHitEE, geo, caloTopology, hotZone_, hotZone_, -500.0, 500.0, (verbose_>1));
231  double e1x1 = spr::eECALmatrix(hotCell, caloHitEB, caloHitEE, geo, caloTopology, 0, 0, -100.0, -100.0,-500.0, 500.0, (verbose_>2));
232  double e3x3 = spr::eECALmatrix(hotCell, caloHitEB, caloHitEE, geo, caloTopology, 1, 1, -100.0, -100.0,-500.0, 500.0, (verbose_>2));
233  double e5x5 = spr::eECALmatrix(hotCell, caloHitEB, caloHitEE, geo, caloTopology, 2, 2, -100.0, -100.0,-500.0, 500.0, (verbose_>2));
234  double p = simTrkItr->momentum().P();
235  double eta = std::abs(simTrkItr->momentum().eta());
236  int etaBin=-1, momBin=-1;
237  for (int ieta=0; ieta<NEtaBins_; ieta++) {
238  if (eta>etaBins_[ieta] && eta<etaBins_[ieta+1] ) etaBin = ieta+1;
239  }
240  for (int ipt=0; ipt<NPBins_; ipt++) {
241  if (p>pBins_[ipt] && p<pBins_[ipt+1] ) momBin = ipt+1;
242  }
243  double r1=-1, r2=-1, r3=-1;
244  if (e3x3 > 0) r1 = e1x1/e3x3;
245  if (e5x5 > 0) {r2 = e1x1/e5x5; r3 = e3x3/e5x5;}
246  if (verbose_>0) {
247  edm::LogVerbatim("IsoTrack") << "ElectronStudy: p " << p << " ["
248  << momBin << "] eta " << eta << " ["
249  << etaBin << "] Cell 0x" << std::hex
250  << isoCell() << std::dec;
251  if (isoCell.subdetId() == EcalBarrel) {
252  edm::LogVerbatim("IsoTrack") << EBDetId(isoCell);
253  } else if (isoCell.subdetId() == EcalEndcap) {
254  edm::LogVerbatim("IsoTrack") << EEDetId(isoCell);
255  }
256  edm::LogVerbatim("IsoTrack") << " e1x1 " << e1x1 << "|" << r1
257  << "|" << r2 << " e3x3 " << e3x3
258  << "|" << r3 << " e5x5 " << e5x5;
259  }
260  if (histos_) {
261  histoR1[0][0]->Fill(r1);
262  histoR2[0][0]->Fill(r2);
263  histoR3[0][0]->Fill(r3);
264  histoE1x1[0][0]->Fill(e1x1/p);
265  histoE3x3[0][0]->Fill(e3x3/p);
266  histoE5x5[0][0]->Fill(e5x5/p);
267  if (momBin>0) {
268  histoR1[momBin][0]->Fill(r1);
269  histoR2[momBin][0]->Fill(r2);
270  histoR3[momBin][0]->Fill(r3);
271  histoE1x1[momBin][0]->Fill(e1x1/p);
272  histoE3x3[momBin][0]->Fill(e3x3/p);
273  histoE5x5[momBin][0]->Fill(e5x5/p);
274  }
275  if (etaBin>0) {
276  histoR1[0][etaBin]->Fill(r1);
277  histoR2[0][etaBin]->Fill(r2);
278  histoR3[0][etaBin]->Fill(r3);
279  histoE1x1[0][etaBin]->Fill(e1x1/p);
280  histoE3x3[0][etaBin]->Fill(e3x3/p);
281  histoE5x5[0][etaBin]->Fill(e5x5/p);
282  if (momBin>0) {
283  histoR1[momBin][etaBin]->Fill(r1);
284  histoR2[momBin][etaBin]->Fill(r2);
285  histoR3[momBin][etaBin]->Fill(r3);
286  histoE1x1[momBin][etaBin]->Fill(e1x1/p);
287  histoE3x3[momBin][etaBin]->Fill(e3x3/p);
288  histoE5x5[momBin][etaBin]->Fill(e5x5/p);
289  }
290  }
291  }
292  }
293  }
294  }
295  }
296 
297 }
298 
299 //define this as a plug-in
RunNumber_t run() const
Definition: EventID.h:39
ElectronStudy(const edm::ParameterSet &ps)
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
std::vector< spr::propagatedTrackID > propagateCALO(edm::Handle< reco::TrackCollection > &trkCollection, const CaloGeometry *geo, const MagneticField *bField, const std::string &theTrackQuality, bool debug=false)
edm::EDGetTokenT< edm::SimTrackContainer > tok_simTk_
ParameterDescriptionBase * addUntracked(U const &iLabel, T const &value)
~ElectronStudy() override
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
TH1F * histoR3[NPBins_+1][NEtaBins_+1]
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< edm::SimVertexContainer > tok_simVtx_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
static const int NEtaBins_
TH1F * histoE3x3[NPBins_+1][NEtaBins_+1]
TH1F * histoR2[NPBins_+1][NEtaBins_+1]
void beginJob() override
int iEvent
Definition: GenABIO.cc:230
bool isAvailable() const
Definition: Service.h:46
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double etaBins_[NEtaBins_+1]
TH1F * histoE1x1[NPBins_+1][NEtaBins_+1]
bool isValid() const
Definition: HandleBase.h:74
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:38
edm::EDGetTokenT< edm::PCaloHitContainer > tok_EBhit_
Definition: DetId.h:18
TH1F * histoE5x5[NPBins_+1][NEtaBins_+1]
const T & get() const
Definition: EventSetup.h:59
DetId hotCrystal(const DetId &detId, edm::Handle< T > &hitsEB, edm::Handle< T > &hitsEE, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, double tMin=-500, double tMax=500, bool debug=false)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< edm::PCaloHitContainer > tok_EEhit_
edm::EventID id() const
Definition: EventBase.h:60
void analyze(edm::Event const &, edm::EventSetup const &) override
static const int NPBins_
T const * product() const
Definition: ESHandle.h:86
double eECALmatrix(const DetId &detId, edm::Handle< T > &hitsEB, edm::Handle< T > &hitsEE, const CaloGeometry *geo, const CaloTopology *caloTopology, int ieta, int iphi, double ebThr=-100, double eeThr=-100, double tMin=-500, double tMax=500, bool debug=false)
TH1F * histoR1[NPBins_+1][NEtaBins_+1]
double pBins_[NPBins_+1]