CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SimHitsValidationHcal.cc
Go to the documentation of this file.
3 
5 
6  g4Label = ps.getUntrackedParameter<std::string>("moduleLabel","g4SimHits");
7  hcalHits = ps.getUntrackedParameter<std::string>("HitCollection","HcalHits");
8  verbose_ = ps.getUntrackedParameter<bool>("Verbose", false);
9 
10  tok_hits_ = consumes<edm::PCaloHitContainer>(edm::InputTag(g4Label,hcalHits));
11 
12  edm::LogInfo("HitsValidationHcal") << "Module Label: " << g4Label << " Hits: "
13  << hcalHits;
14 
15 }
16 
18 
20 
21  edm::LogInfo("HitsValidationHcal") << "Histograms booked";
22  ib.setCurrentFolder("HcalHitsV/SimHitsValidationHcal");
23 
24  //Histograms for Hits
25 
26  std::string divisions[nType]={"HB0","HB1","HE0+z","HE1+z","HE2+z","HE0-z","HE1-z",
27  "HE2-z","HO0","HFL0+z","HFS0+z","HFL1+z","HFS1+z",
28  "HFL2+z","HFS2+z","HFL3+z","HFS3+z","HFL0-z","HFS0-z",
29  "HFL1-z","HFS1-z","HFL2-z","HFS2-z","HFL3-z","HFS3-z"};
30 
31  std::string divisions1[nType]={"HB depth1","HB depth2 ","HE+z depth1","HE+z depth2","HE+z depth3","HE-z depth1","HE-z depth2",
32  "HE-z depth3","HO depth1","HFL+z depth1","HFS+z depth1","HFL+z depth2","HFS+z depth2",
33  "HFL+z depth3","HFS+z depth3","HFL+z depth4","HFS+z depth4","HFL-z depth1","HFS-z depth1",
34  "HFL-z depth2","HFS-z depth2","HFL-z depth3","HFS-z depth3 ","HFL-z depth4","HFS-z depth4"};
35 
36  double etaLow[nType]={-16,-16,16,16,16,-30,-30,-30,-15,29,29,29,29,29,29,29,29,
37  -41,-41,-41,-41,-41,-41,-41,-41};
38  double etaHigh[nType]={16,16,30,30,30,-16,-16,-16,15,41,41,41,41,41,41,41,41,
39  -29,-29,-29,-29,-29,-29,-29,-29};
40  int etaBins[nType]={32,32,14,14,14,14,14,14,30,12,12,12,12,12,12,12,12,
41  12,12,12,12,12,12,12,12};
42  char name[40], title[100];
43 
44  for (int i=0; i<nType; ++i) {
45  sprintf (name, "HcalHitEta%s", divisions[i].c_str());
46  sprintf (title, "Hit energy as a function of eta tower index in %s", divisions1[i].c_str());
47  meHcalHitEta_[i] = ib.book1D(name, title, etaBins[i], etaLow[i], etaHigh[i]);
48 
49  sprintf (name, "HcalHitTimeAEta%s", divisions[i].c_str());
50  sprintf (title, "Hit time as a function of eta tower index in %s", divisions1[i].c_str());
51  meHcalHitTimeEta_[i] = ib.book1D(name, title, etaBins[i], etaLow[i], etaHigh[i]);
52 
53  sprintf (name, "HcalHitE25%s", divisions[i].c_str());
54  sprintf (title, "Energy in time window 0 to 25 for a tower in %s", divisions1[i].c_str());
55  meHcalEnergyl25_[i] = ib.book2D(name, title, etaBins[i], etaLow[i], etaHigh[i], 72, 0., 72.);
56 
57  sprintf (name, "HcalHitE50%s", divisions[i].c_str());
58  sprintf (title, "Energy in time window 0 to 50 for a tower in %s", divisions1[i].c_str());
59  meHcalEnergyl50_[i] = ib.book2D(name, title, etaBins[i], etaLow[i], etaHigh[i], 72, 0., 72.);
60 
61  sprintf (name, "HalHitE100%s", divisions[i].c_str());
62  sprintf (title, "Energy in time window 0 to 100 for a tower in %s", divisions1[i].c_str());
63  meHcalEnergyl100_[i] = ib.book2D(name, title, etaBins[i], etaLow[i], etaHigh[i], 72, 0., 72.);
64 
65  sprintf (name, "HcalHitE250%s", divisions[i].c_str());
66  sprintf (title, "Energy in time window 0 to 250 for a tower in %s", divisions1[i].c_str());
67  meHcalEnergyl250_[i] = ib.book2D(name, title, etaBins[i], etaLow[i], etaHigh[i], 72, 0., 72.);
68  }
69 
70  sprintf (name, "Energy_HB");
71  meEnergy_HB = ib.book1D(name, name, 100,0,1);
72  sprintf (name, "Energy_HE");
73  meEnergy_HE = ib.book1D(name, name, 100,0,1);
74  sprintf (name, "Energy_HO");
75  meEnergy_HO = ib.book1D(name, name, 100,0,1);
76  sprintf (name, "Energy_HF");
77  meEnergy_HF = ib.book1D(name, name, 100,0,50);
78 
79  sprintf (name, "Time_HB");
80  metime_HB = ib.book1D(name, name, 300,-150,150);
81  sprintf (name, "Time_HE");
82  metime_HE = ib.book1D(name, name, 300,-150,150);
83  sprintf (name, "Time_HO");
84  metime_HO = ib.book1D(name, name, 300,-150, 150);
85  sprintf (name, "Time_HF");
86  metime_HF = ib.book1D(name, name, 300,-150,150);
87 
88  sprintf (name, "Time_Enweighted_HB");
89  metime_enweighted_HB = ib.book1D(name, name, 300,-150,150);
90  sprintf (name, "Time_Enweighted_HE");
91  metime_enweighted_HE = ib.book1D(name, name, 300,-150,150);
92  sprintf (name, "Time_Enweighted_HO");
93  metime_enweighted_HO = ib.book1D(name, name, 300,-150, 150);
94  sprintf (name, "Time_Enweighted_HF");
95  metime_enweighted_HF = ib.book1D(name, name, 300,-150,150);
96 }
97 
99 
100 
101  LogDebug("SimHitsValidationHcal") << "Run = " << e.id().run() << " Event = "
102  << e.id().event();
103 
104  std::vector<PCaloHit> caloHits;
106 
107  bool getHits = false;
108  e.getByToken(tok_hits_,hitsHcal);
109  if (hitsHcal.isValid()) getHits = true;
110 
111  LogDebug("SimHitsValidationHcal") << "SimHitsValidationHcal.: Input flags Hits " << getHits;
112 
113  if (getHits) {
114  caloHits.insert(caloHits.end(),hitsHcal->begin(),hitsHcal->end());
115  LogDebug("SimHitsValidationHcal") << "SimHitsValidationHcal: Hit buffer "
116  << caloHits.size();
117  analyzeHits (caloHits);
118  }
119 }
120 
121 void SimHitsValidationHcal::analyzeHits (std::vector<PCaloHit>& hits) {
122 
123  int nHit = hits.size();
124  double entotHB = 0, entotHE = 0, entotHF = 0, entotHO = 0;
125  double timetotHB = 0, timetotHE = 0, timetotHF = 0, timetotHO = 0;
126  int nHB=0, nHE=0, nHO=0, nHF=0;
127 
128  std::map<std::pair<HcalDetId,int>,energysum> map_try;
129  map_try.clear();
130 
131  std::map<std::pair<HcalDetId,int>,energysum>::iterator itr;
132 
133  for (int i=0; i<nHit; i++) {
134  double energy = hits[i].energy();
135  double time = hits[i].time();
136  unsigned int id_ = hits[i].id();
137  HcalDetId id = HcalDetId(id_);
138  int itime = (int)(time);
139  int det = id.det();
140  int subdet = id.subdet();
141  int depth = id.depth();
142  int eta = id.ieta();
143  int phi = id.iphi();
144  unsigned int dep = hits[i].depth();
145 
146  int type =-1;
147  if (subdet == static_cast<int>(HcalBarrel)) {
148  entotHB += energy;
149  timetotHB += time;
150  nHB++;
151  type = depth-1;
152  } else if (subdet == static_cast<int>(HcalEndcap)) {
153  entotHE += energy;
154  timetotHE += time;
155  nHE++;
156  type = depth+1;
157  if (eta < 0) type += 3;
158  } else if (subdet == static_cast<int>(HcalOuter)) {
159  entotHO += energy;
160  timetotHO += time;
161  nHO++;
162  type = 8;
163  } else if (subdet == static_cast<int>(HcalForward)) {
164  entotHF += energy;
165  timetotHF += time;
166  nHF++;
167  type = depth+8+2*dep;
168  if (eta < 0) type += 8;
169  }
170 
171  std::pair<HcalDetId,int> id0(id,type);
172  // LogDebug("HitsValidationHcal") << "Id and type " << id << ":" << type;
173  energysum ensum;
174  if (map_try.count(id0) != 0) ensum = map_try[id0];
175  if (itime<250) {
176  ensum.e250 += energy;
177  if (itime<100) {
178  ensum.e100 += energy;
179  if (itime<50) {
180  ensum.e50 += energy;
181  if (itime<25) ensum.e25 += energy;
182  }
183  }
184  }
185  map_try[id0] = ensum;
186 
187  LogDebug("SimHitsValidationHcal") << "Hit[" << i << "] ID " << std::hex << id_
188  << std::dec << " " << id << std::dec<< " Det " << det << " Sub "
189  << subdet << " depth " << depth << " depthX "
190  << dep << " Eta " << eta << " Phi " << phi
191  << " E " << energy << " time " << time
192  << " type " << type;
193 
194  // LogDebug("HitsValidationHcal") << "Hit[" << i << "] ID " << std::hex << id_ << "ID="<<std::dec << id << " Det " << det << " Sub " << subdet << " depth " << depth << " depthX " << dep << " Eta " << eta << " Phi " << phi << " E " << energy << " time " << time <<" itime " << itime << " type " << type;
195 
196  double etax = eta - 0.5;
197  if (eta < 0) etax += 1;
198  if ( type >= 0) {
199  meHcalHitEta_[type]->Fill(etax,energy);
200  meHcalHitTimeEta_[type]->Fill(etax,time);
201  }
202  }
203 
204  meEnergy_HB->Fill(entotHB);
205  meEnergy_HE->Fill(entotHE);
206  meEnergy_HF->Fill(entotHF);
207  meEnergy_HO->Fill(entotHO);
208 
209  metime_HB->Fill(timetotHB);
210  metime_HE->Fill(timetotHE);
211  metime_HF->Fill(timetotHF);
212  metime_HO->Fill(timetotHO);
213 
214  metime_enweighted_HB->Fill(timetotHB,entotHB);
215  metime_enweighted_HE->Fill(timetotHE,entotHE);
216  metime_enweighted_HF->Fill(timetotHF,entotHF);
217  metime_enweighted_HO->Fill(timetotHO,entotHO);
218 
219 
220  for ( itr = map_try.begin() ; itr != map_try.end(); ++itr) {
221  if ( (*itr).first.second >= 0) {
222  HcalDetId id= (*itr).first.first;
223  energysum ensum= (*itr).second;
224  int eta = id.ieta();
225  int phi = id.iphi();
226  double etax= eta-0.5;
227  double phix= phi-0.5;
228  meHcalEnergyl25_[(*itr).first.second]->Fill(etax,phix,ensum.e25);
229  meHcalEnergyl50_[(*itr).first.second]->Fill(etax,phix,ensum.e50);
230  meHcalEnergyl100_[(*itr).first.second]->Fill(etax,phix,ensum.e100);
231  meHcalEnergyl250_[(*itr).first.second]->Fill(etax,phix,ensum.e250);
232 
233  //LogDebug("HitsValidationHcal") <<" energy of tower ="<<(*itr).first<<"in time 25ns is== "<<(*itr).second.e25<<"in time 25-50ns=="<<(*itr).second.e50<<"in time 50-100ns=="<<(*itr).second.e100<<"in time 100-250 ns=="<<(*itr).second.e250;
234  }
235  }
236 
237 }
238 
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:39
SimHitsValidationHcal(const edm::ParameterSet &ps)
type
Definition: HCALResponse.h:21
MonitorElement * meEnergy_HF
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
int ib
Definition: cuy.py:660
MonitorElement * meHcalHitEta_[nType]
MonitorElement * meHcalEnergyl50_[nType]
MonitorElement * metime_enweighted_HB
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:446
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &)
void analyze(const edm::Event &e, const edm::EventSetup &c)
T eta() const
MonitorElement * meHcalHitTimeEta_[nType]
void Fill(long long x)
MonitorElement * meEnergy_HB
MonitorElement * meHcalEnergyl250_[nType]
MonitorElement * metime_enweighted_HF
MonitorElement * meHcalEnergyl25_[nType]
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:115
bool isValid() const
Definition: HandleBase.h:76
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:133
MonitorElement * meEnergy_HE
edm::EventID id() const
Definition: EventBase.h:56
MonitorElement * meHcalEnergyl100_[nType]
MonitorElement * metime_enweighted_HE
MonitorElement * metime_enweighted_HO
MonitorElement * meEnergy_HO
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hits_
Definition: Run.h:41
Definition: DDAxes.h:10
void analyzeHits(std::vector< PCaloHit > &)