CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalSimHitStudy.cc
Go to the documentation of this file.
3 
5 
7 
8  g4Label = ps.getUntrackedParameter<std::string>("moduleLabel","g4SimHits");
9  hcalHits = ps.getUntrackedParameter<std::string>("HitCollection","HcalHits");
10  outFile_ = ps.getUntrackedParameter<std::string>("outputFile", "hcHit.root");
11  verbose_ = ps.getUntrackedParameter<bool>("Verbose", false);
12  checkHit_= true;
13 
14  edm::LogInfo("HcalSim") << "Module Label: " << g4Label << " Hits: "
15  << hcalHits << " / "<< checkHit_
16  << " Output: " << outFile_;
17 
19  if (dbe_) {
20  if (verbose_) {
21  dbe_->setVerbose(1);
22  sleep (3);
24  } else {
25  dbe_->setVerbose(0);
26  }
27  }
28 }
29 
31 
33 
34  if (dbe_) {
35  dbe_->setCurrentFolder("HcalHitsV/HcalSimHitsTask");
36 
37  //Histograms for Hits
38  if (checkHit_) {
39  meAllNHit_ = dbe_->book1D("Hit01","Number of Hits in HCal",1000,0.,1000.);
40  meBadDetHit_= dbe_->book1D("Hit02","Hits with wrong Det", 100,0.,100.);
41  meBadSubHit_= dbe_->book1D("Hit03","Hits with wrong Subdet",100,0.,100.);
42  meBadIdHit_ = dbe_->book1D("Hit04","Hits with wrong ID", 100,0.,100.);
43  meHBNHit_ = dbe_->book1D("Hit05","Number of Hits in HB",1000,0.,1000.);
44  meHENHit_ = dbe_->book1D("Hit06","Number of Hits in HE",1000,0.,1000.);
45  meHONHit_ = dbe_->book1D("Hit07","Number of Hits in HO",1000,0.,1000.);
46  meHFNHit_ = dbe_->book1D("Hit08","Number of Hits in HF",1000,0.,1000.);
47  meDetectHit_= dbe_->book1D("Hit09","Detector ID", 50,0.,50.);
48  meSubdetHit_= dbe_->book1D("Hit10","Subdetectors in HCal", 50,0.,50.);
49  meDepthHit_ = dbe_->book1D("Hit11","Depths in HCal", 20,0.,20.);
50  meEtaHit_ = dbe_->book1D("Hit12","Eta in HCal", 100,-50.,50.);
51  mePhiHit_ = dbe_->book1D("Hit13","Phi in HCal", 100,0.,100.);
52  meEnergyHit_= dbe_->book1D("Hit14","Energy in HCal", 100,0.,1.);
53  meTimeHit_ = dbe_->book1D("Hit15","Time in HCal", 100,0.,400.);
54  meTimeWHit_ = dbe_->book1D("Hit16","Time in HCal (E wtd)", 100,0.,400.);
55  meHBDepHit_ = dbe_->book1D("Hit17","Depths in HB", 20,0.,20.);
56  meHEDepHit_ = dbe_->book1D("Hit18","Depths in HE", 20,0.,20.);
57  meHODepHit_ = dbe_->book1D("Hit19","Depths in HO", 20,0.,20.);
58  meHFDepHit_ = dbe_->book1D("Hit20","Depths in HF", 20,0.,20.);
59  meHBEtaHit_ = dbe_->book1D("Hit21","Eta in HB", 100,-50.,50.);
60  meHEEtaHit_ = dbe_->book1D("Hit22","Eta in HE", 100,-50.,50.);
61  meHOEtaHit_ = dbe_->book1D("Hit23","Eta in HO", 100,-50.,50.);
62  meHFEtaHit_ = dbe_->book1D("Hit24","Eta in HF", 100,-50.,50.);
63  meHBPhiHit_ = dbe_->book1D("Hit25","Phi in HB", 100,0.,100.);
64  meHEPhiHit_ = dbe_->book1D("Hit26","Phi in HE", 100,0.,100.);
65  meHOPhiHit_ = dbe_->book1D("Hit27","Phi in HO", 100,0.,100.);
66  meHFPhiHit_ = dbe_->book1D("Hit28","Phi in HF", 100,0.,100.);
67  meHBEneHit_ = dbe_->book1D("Hit29","Energy in HB", 100,0.,1.);
68  meHEEneHit_ = dbe_->book1D("Hit30","Energy in HE", 100,0.,1.);
69  meHOEneHit_ = dbe_->book1D("Hit31","Energy in HO", 100,0.,1.);
70  meHFEneHit_ = dbe_->book1D("Hit32","Energy in HF", 100,0.,100.);
71  meHBTimHit_ = dbe_->book1D("Hit33","Time in HB", 100,0.,400.);
72  meHETimHit_ = dbe_->book1D("Hit34","Time in HE", 100,0.,400.);
73  meHOTimHit_ = dbe_->book1D("Hit35","Time in HO", 100,0.,400.);
74  meHFTimHit_ = dbe_->book1D("Hit36","Time in HF", 100,0.,400.);
75  meHBEneHit2_ = dbe_->book1D("Hit37","Energy in HB 2", 100,0.,0.0001);
76  meHEEneHit2_ = dbe_->book1D("Hit38","Energy in HE 2", 100,0.,0.0001);
77  meHOEneHit2_ = dbe_->book1D("Hit39","Energy in HO 2", 100,0.,0.0001);
78  meHFEneHit2_ = dbe_->book1D("Hit40","Energy in HF 2", 100,0.,0.0001);
79  meHBL10Ene_ = dbe_->book1D("Hit41","Log10Energy in HB", 140, -10., 4. );
80  meHEL10Ene_ = dbe_->book1D("Hit42","Log10Energy in HE", 140, -10., 4. );
81  meHFL10Ene_ = dbe_->book1D("Hit43","Log10Energy in HF", 140, -10., 4. );
82  meHOL10Ene_ = dbe_->book1D("Hit44","Log10Energy in HO", 140, -10., 4. );
83  meHBL10EneP_ = dbe_->bookProfile("Hit45","Log10Energy in HB vs Hit contribution", 140, -10., 4., 100, 0., 1. );
84  meHEL10EneP_ = dbe_->bookProfile("Hit46","Log10Energy in HE vs Hit contribution", 140, -10., 4., 100, 0., 1. );
85  meHFL10EneP_ = dbe_->bookProfile("Hit47","Log10Energy in HF vs Hit contribution", 140, -10., 4., 100, 0., 1. );
86  meHOL10EneP_ = dbe_->bookProfile("Hit48","Log10Energy in HO vs Hit contribution", 140, -10., 4., 100, 0., 1. );
87  }
88  }
89 }
90 
92  if (dbe_ && outFile_.size() > 0) dbe_->save(outFile_);
93 }
94 
96 
97  LogDebug("HcalSim") << "Run = " << e.id().run() << " Event = "
98  << e.id().event();
99 
100  std::vector<PCaloHit> caloHits;
102 
103  bool getHits = false;
104  if (checkHit_) {
105  e.getByLabel(g4Label,hcalHits,hitsHcal);
106  if (hitsHcal.isValid()) getHits = true;
107  }
108 
109  LogDebug("HcalSim") << "HcalValidation: Input flags Hits " << getHits;
110 
111  if (getHits) {
112  caloHits.insert(caloHits.end(),hitsHcal->begin(),hitsHcal->end());
113  LogDebug("HcalSim") << "HcalValidation: Hit buffer "
114  << caloHits.size();
115  analyzeHits (caloHits);
116  }
117 }
118 
119 void HcalSimHitStudy::analyzeHits (std::vector<PCaloHit>& hits) {
120 
121  int nHit = hits.size();
122  int nHB=0, nHE=0, nHO=0, nHF=0, nBad1=0, nBad2=0, nBad=0;
123  std::vector<double> encontHB(140, 0.);
124  std::vector<double> encontHE(140, 0.);
125  std::vector<double> encontHF(140, 0.);
126  std::vector<double> encontHO(140, 0.);
127  double entotHB = 0, entotHE = 0, entotHF = 0, entotHO = 0;
128 
129  for (int i=0; i<nHit; i++) {
130  double energy = hits[i].energy();
131  double log10en = log10(energy);
132  int log10i = int( (log10en+10.)*10. );
133  double time = hits[i].time();
134  unsigned int id_ = hits[i].id();
135  HcalDetId id = HcalDetId(id_);
136  int det = id.det();
137  int subdet = id.subdet();
138  int depth = id.depth();
139  int eta = id.ieta();
140  int phi = id.iphi();
141  LogDebug("HcalSim") << "Hit[" << i << "] ID " << std::hex << id_
142  << std::dec << " Det " << det << " Sub "
143  << subdet << " depth " << depth << " Eta " << eta
144  << " Phi " << phi << " E " << energy << " time "
145  << time;
146  if (det == 4) { // Check DetId.h
147  if (subdet == static_cast<int>(HcalBarrel)) nHB++;
148  else if (subdet == static_cast<int>(HcalEndcap)) nHE++;
149  else if (subdet == static_cast<int>(HcalOuter)) nHO++;
150  else if (subdet == static_cast<int>(HcalForward)) nHF++;
151  else { nBad++; nBad2++;}
152  } else { nBad++; nBad1++;}
153  if (dbe_) {
154  meDetectHit_->Fill(double(det));
155  if (det == 4) {
156  meSubdetHit_->Fill(double(subdet));
157  meDepthHit_->Fill(double(depth));
158  meEtaHit_->Fill(double(eta));
159  mePhiHit_->Fill(double(phi));
160  meEnergyHit_->Fill(energy);
161  meTimeHit_->Fill(time);
162  meTimeWHit_->Fill(double(time),energy);
163  if (subdet == static_cast<int>(HcalBarrel)) {
164  meHBDepHit_->Fill(double(depth));
165  meHBEtaHit_->Fill(double(eta));
166  meHBPhiHit_->Fill(double(phi));
167  meHBEneHit_->Fill(energy);
168  meHBEneHit2_->Fill(energy);
169  meHBTimHit_->Fill(time);
170  meHBL10Ene_->Fill(log10en);
171  if( log10i >=0 && log10i < 140 ) encontHB[log10i] += energy;
172  entotHB += energy;
173  } else if (subdet == static_cast<int>(HcalEndcap)) {
174  meHEDepHit_->Fill(double(depth));
175  meHEEtaHit_->Fill(double(eta));
176  meHEPhiHit_->Fill(double(phi));
177  meHEEneHit_->Fill(energy);
178  meHEEneHit2_->Fill(energy);
179  meHETimHit_->Fill(time);
180  meHEL10Ene_->Fill(log10en);
181  if( log10i >=0 && log10i < 140 ) encontHE[log10i] += energy;
182  entotHE += energy;
183  } else if (subdet == static_cast<int>(HcalOuter)) {
184  meHODepHit_->Fill(double(depth));
185  meHOEtaHit_->Fill(double(eta));
186  meHOPhiHit_->Fill(double(phi));
187  meHOEneHit_->Fill(energy);
188  meHOEneHit2_->Fill(energy);
189  meHOTimHit_->Fill(time);
190  meHOL10Ene_->Fill(log10en);
191  if( log10i >=0 && log10i < 140 ) encontHO[log10i] += energy;
192  entotHO += energy;
193  } else if (subdet == static_cast<int>(HcalForward)) {
194  meHFDepHit_->Fill(double(depth));
195  meHFEtaHit_->Fill(double(eta));
196  meHFPhiHit_->Fill(double(phi));
197  meHFEneHit_->Fill(energy);
198  meHFEneHit2_->Fill(energy);
199  meHFTimHit_->Fill(time);
200  meHFL10Ene_->Fill(log10en);
201  if( log10i >=0 && log10i < 140 ) encontHF[log10i] += energy;
202  entotHF += energy;
203  }
204  }
205  }
206  }
207  if( entotHB != 0 ) for( int i=0; i<140; i++ ) meHBL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHB[i]/entotHB );
208  if( entotHE != 0 ) for( int i=0; i<140; i++ ) meHEL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHE[i]/entotHE );
209  if( entotHF != 0 ) for( int i=0; i<140; i++ ) meHFL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHF[i]/entotHF );
210  if( entotHO != 0 ) for( int i=0; i<140; i++ ) meHOL10EneP_->Fill( -10.+(float(i)+0.5)/10., encontHO[i]/entotHO );
211 
212  if (dbe_) {
213  meAllNHit_->Fill(double(nHit));
214  meBadDetHit_->Fill(double(nBad1));
215  meBadSubHit_->Fill(double(nBad2));
216  meBadIdHit_->Fill(double(nBad));
217  meHBNHit_->Fill(double(nHB));
218  meHENHit_->Fill(double(nHE));
219  meHONHit_->Fill(double(nHO));
220  meHFNHit_->Fill(double(nHF));
221  }
222  LogDebug("HcalSim") << "HcalSimHitStudy::analyzeHits: HB " << nHB
223  << " HE " << nHE << " HO " << nHO << " HF " << nHF
224  << " Bad " << nBad << " All " << nHit;
225 
226 }
227 
#define LogDebug(id)
RunNumber_t run() const
Definition: EventID.h:42
MonitorElement * mePhiHit_
EventNumber_t event() const
Definition: EventID.h:44
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
MonitorElement * meHEDepHit_
MonitorElement * meHOL10Ene_
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:722
MonitorElement * meHOPhiHit_
MonitorElement * meDetectHit_
MonitorElement * meHFNHit_
MonitorElement * meHFEtaHit_
MonitorElement * meHFL10EneP_
MonitorElement * meHFEneHit2_
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2118
MonitorElement * meHOL10EneP_
MonitorElement * meHBNHit_
MonitorElement * meHFTimHit_
T eta() const
void sleep(Duration_t)
Definition: Utils.h:163
MonitorElement * meHBEneHit2_
MonitorElement * meHONHit_
void Fill(long long x)
MonitorElement * meAllNHit_
MonitorElement * meHBL10Ene_
MonitorElement * meBadDetHit_
MonitorElement * meHFL10Ene_
MonitorElement * meHETimHit_
MonitorElement * meTimeHit_
void analyzeHits(std::vector< PCaloHit > &)
MonitorElement * meHFDepHit_
MonitorElement * meHBTimHit_
MonitorElement * meTimeWHit_
MonitorElement * bookProfile(const char *name, const char *title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const char *option="s")
Definition: DQMStore.cc:1036
void setVerbose(unsigned level)
Definition: DQMStore.cc:398
bool isValid() const
Definition: HandleBase.h:76
MonitorElement * meHBDepHit_
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
void analyze(const edm::Event &e, const edm::EventSetup &c)
MonitorElement * meBadSubHit_
MonitorElement * meHENHit_
MonitorElement * meBadIdHit_
MonitorElement * meHEL10EneP_
MonitorElement * meEnergyHit_
std::string hcalHits
MonitorElement * meHOEtaHit_
MonitorElement * meEtaHit_
MonitorElement * meHBEneHit_
MonitorElement * meHOTimHit_
MonitorElement * meHBL10EneP_
edm::EventID id() const
Definition: EventBase.h:56
MonitorElement * meHEEneHit2_
MonitorElement * meHFPhiHit_
std::string outFile_
MonitorElement * meDepthHit_
MonitorElement * meHBEtaHit_
MonitorElement * meHODepHit_
MonitorElement * meHEL10Ene_
void showDirStructure(void) const
Definition: DQMStore.cc:2766
MonitorElement * meHEEneHit_
std::string g4Label
MonitorElement * meHEPhiHit_
HcalSimHitStudy(const edm::ParameterSet &ps)
MonitorElement * meSubdetHit_
MonitorElement * meHBPhiHit_
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:434
MonitorElement * meHFEneHit_
MonitorElement * meHEEtaHit_
MonitorElement * meHOEneHit_
MonitorElement * meHOEneHit2_
Definition: DDAxes.h:10