CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalSimHitsValidation.cc
Go to the documentation of this file.
3 
4 
6  // DQM ROOT output
7  outputFile_ = conf.getUntrackedParameter<std::string>("outputFile", "myfile.root");
8 
9  // register for data access
10  tok_evt_ = consumes<edm::HepMCProduct>(edm::InputTag("generator"));
11  tok_hcal_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits","HcalHits"));
12  tok_ecalEB_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits","EcalHitsEB"));
13  tok_ecalEE_ = consumes<edm::PCaloHitContainer>(edm::InputTag("g4SimHits","EcalHitsEE"));
14 
15  if ( outputFile_.size() != 0 ) { edm::LogInfo("OutputInfo") << " Hcal SimHit Task histograms will be saved to '" << outputFile_.c_str() << "'";
16  } else {
17  edm::LogInfo("OutputInfo") << " Hcal SimHit Task histograms will NOT be saved";
18  }
19 
20  nevtot = 0;
21 
22  dbe_ = 0;
23  // get hold of back-end interface
25 
26  Char_t histo[200];
27 
28  if ( dbe_ ) {
29  dbe_->setCurrentFolder("HcalSimHitsV/HcalSimHitTask");
30 
31  // General counters
32  sprintf (histo, "N_HB" );
33  Nhb = dbe_->book1D(histo, histo, 2600,0.,2600.);
34  sprintf (histo, "N_HE" );
35  Nhe = dbe_->book1D(histo, histo, 2600,0.,2600.);
36  sprintf (histo, "N_HO" );
37  Nho = dbe_->book1D(histo, histo, 2200,0.,2200.);
38  sprintf (histo, "N_HF" );
39  Nhf = dbe_->book1D(histo, histo, 1800,0., 1800.);
40 
41  //Mean energy vs iEta TProfiles
42  sprintf (histo, "emean_vs_ieta_HB1" );
43  emean_vs_ieta_HB1 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s");
44  sprintf (histo, "emean_vs_ieta_HB2" );
45  emean_vs_ieta_HB2 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s");
46  sprintf (histo, "emean_vs_ieta_HE1" );
47  emean_vs_ieta_HE1 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10. ,2000., "s" );
48  sprintf (histo, "emean_vs_ieta_HE2" );
49  emean_vs_ieta_HE2 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s");
50  sprintf (histo, "emean_vs_ieta_HE3" );
51  emean_vs_ieta_HE3 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
52  sprintf (histo, "emean_vs_ieta_HO" );
53  emean_vs_ieta_HO = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
54  sprintf (histo, "emean_vs_ieta_HF1" );
55  emean_vs_ieta_HF1 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
56  sprintf (histo, "emean_vs_ieta_HF2" );
57  emean_vs_ieta_HF2 = dbe_->bookProfile(histo, histo, 82, -41., 41., 2010, -10., 2000., "s" );
58 
59  //Occupancy vs. iEta TH1Fs
60  sprintf (histo, "occupancy_vs_ieta_HB1" );
61  occupancy_vs_ieta_HB1 = dbe_->book1D(histo, histo, 82, -41., 41.);
62  sprintf (histo, "occupancy_vs_ieta_HB2" );
63  occupancy_vs_ieta_HB2 = dbe_->book1D(histo, histo, 82, -41., 41.);
64  sprintf (histo, "occupancy_vs_ieta_HE1" );
65  occupancy_vs_ieta_HE1 = dbe_->book1D(histo, histo, 82, -41., 41.);
66  sprintf (histo, "occupancy_vs_ieta_HE2" );
67  occupancy_vs_ieta_HE2 = dbe_->book1D(histo, histo, 82, -41., 41.);
68  sprintf (histo, "occupancy_vs_ieta_HE3" );
69  occupancy_vs_ieta_HE3 = dbe_->book1D(histo, histo, 82, -41., 41.);
70  sprintf (histo, "occupancy_vs_ieta_HO" );
71  occupancy_vs_ieta_HO = dbe_->book1D(histo, histo, 82, -41., 41.);
72  sprintf (histo, "occupancy_vs_ieta_HF1" );
73  occupancy_vs_ieta_HF1 = dbe_->book1D(histo, histo, 82, -41., 41.);
74  sprintf (histo, "occupancy_vs_ieta_HF2" );
75  occupancy_vs_ieta_HF2 = dbe_->book1D(histo, histo, 82, -41., 41.);
76 
77  //Energy spectra
78  sprintf (histo, "HcalSimHitTask_energy_of_simhits_HB" ) ;
79  meSimHitsEnergyHB = dbe_->book1D(histo, histo, 510 , -0.1 , 5.);
80 
81  sprintf (histo, "HcalSimHitTask_energy_of_simhits_HE" ) ;
82  meSimHitsEnergyHE = dbe_->book1D(histo, histo, 510, -0.02, 2.);
83 
84  sprintf (histo, "HcalSimHitTask_energy_of_simhits_HO" ) ;
85  meSimHitsEnergyHO = dbe_->book1D(histo, histo, 510 , -0.1 , 5.);
86 
87  sprintf (histo, "HcalSimHitTask_energy_of_simhits_HF" ) ;
88  meSimHitsEnergyHF = dbe_->book1D(histo, histo, 1010 , -5. , 500.);
89 
90  //Energy in Cone
91  sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths");
92  meEnConeEtaProfile = dbe_->bookProfile(histo, histo, 82, -41., 41., 210, -10., 200.);
93 
94  sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths_E");
95  meEnConeEtaProfile_E = dbe_->bookProfile(histo, histo, 82, -41., 41., 210, -10., 200.);
96 
97  sprintf (histo, "HcalSimHitTask_En_simhits_cone_profile_vs_ieta_all_depths_EH");
98  meEnConeEtaProfile_EH = dbe_->bookProfile(histo, histo, 82, -41., 41., 210, -10., 200.);
99 
100  } //end-of if(_dbe)
101 
102 }
103 
104 
106 
108  //before check that histos are there....
109 
110  // check if ME still there (and not killed by MEtoEDM for memory saving)
111  if( dbe_ )
112  {
113  // check existence of first histo in the list
114  if (! dbe_->get("HcalSimHitsV/HcalSimHitTask/N_HB")) return;
115  }
116  else
117  return;
118 
119  //======================================
120 
121  for (int i = 1; i <= occupancy_vs_ieta_HB1->getNbinsX(); i++){
122 
123  int ieta = i - 42; // -41 -1, 0 40
124  if(ieta >=0 ) ieta +=1; // -41 -1, 1 41 - to make it detector-like
125 
126  float phi_factor;
127 
128  if (fabs(ieta) <= 20) phi_factor = 72.;
129  else if (fabs(ieta) < 40) phi_factor = 36.;
130  else phi_factor = 18.;
131 
132  float cnorm;
133 
134  //Occupancy vs. iEta TH1Fs
135  cnorm = occupancy_vs_ieta_HB1->getBinContent(i) / (phi_factor * nevtot);
137  cnorm = occupancy_vs_ieta_HB2->getBinContent(i) / (phi_factor * nevtot);
139 
140  cnorm = occupancy_vs_ieta_HE1->getBinContent(i) / (phi_factor * nevtot);
142  cnorm = occupancy_vs_ieta_HE2->getBinContent(i) / (phi_factor * nevtot);
144  cnorm = occupancy_vs_ieta_HE3->getBinContent(i) / (phi_factor * nevtot);
146 
147  cnorm = occupancy_vs_ieta_HO->getBinContent(i) / (phi_factor * nevtot);
149 
150  cnorm = occupancy_vs_ieta_HF1->getBinContent(i) / (phi_factor * nevtot);
152  cnorm = occupancy_vs_ieta_HF2->getBinContent(i) / (phi_factor * nevtot);
154  }
155 
156 
157  if ( outputFile_.size() != 0 && dbe_ ) dbe_->save(outputFile_);
158 }
159 
160 
162 
164 
165  using namespace edm;
166  using namespace std;
167 
168  //===========================================================================
169  // Getting SimHits
170  //===========================================================================
171 
172  double phi_MC = -999.; // phi of initial particle from HepMC
173  double eta_MC = -999.; // eta of initial particle from HepMC
174 
176  ev.getByToken(tok_evt_,evtMC); // generator in late 310_preX
177  if (!evtMC.isValid()) {
178  std::cout << "no HepMCProduct found" << std::endl;
179  }
180 
181  // MC particle with highest pt is taken as a direction reference
182  double maxPt = -99999.;
183  int npart = 0;
184 
185  const HepMC::GenEvent * myGenEvent = evtMC->GetEvent();
186  for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
187  p != myGenEvent->particles_end(); ++p ) {
188  double phip = (*p)->momentum().phi();
189  double etap = (*p)->momentum().eta();
190  double pt = (*p)->momentum().perp();
191  if(pt > maxPt) {npart++; maxPt = pt; phi_MC = phip; eta_MC = etap; }
192  }
193 
194  double partR = 0.3;
195 
196 
197  //Hcal SimHits
198 
199  //Approximate calibration constants
200  const float calib_HB = 120.;
201  const float calib_HE = 190.;
202  const float calib_HF1 = 1.0/0.383;
203  const float calib_HF2 = 1.0/0.368;
204 
206  ev.getByToken(tok_hcal_,hcalHits);
207  const PCaloHitContainer * SimHitResult = hcalHits.product () ;
208 
209  float eta_diff;
210  float etaMax = 9999;
211  int ietaMax = 0;
212 
213  double HcalCone = 0;
214 
215  c.get<CaloGeometryRecord>().get (geometry);
216 
217  for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResult->begin () ; SimHits != SimHitResult->end(); ++SimHits) {
218  HcalDetId cell(SimHits->id());
219  const CaloCellGeometry* cellGeometry = geometry->getSubdetectorGeometry (cell)->getGeometry (cell);
220  double etaS = cellGeometry->getPosition().eta () ;
221  double phiS = cellGeometry->getPosition().phi () ;
222  double en = SimHits->energy();
223 
224  int sub = cell.subdet();
225  int depth = cell.depth();
226  double ieta = cell.ieta();
227 
228  //Energy in Cone
229  double r = dR(eta_MC, phi_MC, etaS, phiS);
230 
231  if (r < partR){
232  eta_diff = fabs(eta_MC - etaS);
233  if(eta_diff < etaMax) {
234  etaMax = eta_diff;
235  ietaMax = cell.ieta();
236  }
237  //Approximation of calibration
238  if (sub == 1) HcalCone += en*calib_HB;
239  else if (sub == 2) HcalCone += en*calib_HE;
240  else if (sub == 4 && depth == 1) HcalCone += en*calib_HF1;
241  else if (sub == 4 && depth == 2) HcalCone += en*calib_HF2;
242  }
243 
244  //Account for lack of ieta = 0
245  if (ieta > 0) ieta--;
246 
247  //HB
248  if (sub == 1){
249  meSimHitsEnergyHB->Fill(en);
250  if (depth == 1){
251  emean_vs_ieta_HB1->Fill(double(ieta), en);
252  occupancy_vs_ieta_HB1->Fill(double(ieta));
253  }
254  if (depth == 2){
255  emean_vs_ieta_HB2->Fill(double(ieta), en);
256  occupancy_vs_ieta_HB2->Fill(double(ieta));
257  }
258  }
259  //HE
260  if (sub == 2){
261  meSimHitsEnergyHE->Fill(en);
262  if (depth == 1){
263  emean_vs_ieta_HE1->Fill(double(ieta), en);
264  occupancy_vs_ieta_HE1->Fill(double(ieta));
265  }
266  if (depth == 2){
267  emean_vs_ieta_HE2->Fill(double(ieta), en);
268  occupancy_vs_ieta_HE2->Fill(double(ieta));
269  }
270  if (depth == 3){
271  emean_vs_ieta_HE3->Fill(double(ieta), en);
272  occupancy_vs_ieta_HE3->Fill(double(ieta));
273  }
274  }
275  //HO
276  if (sub == 3){
277  meSimHitsEnergyHO->Fill(en);
278  emean_vs_ieta_HO->Fill(double(ieta), en);
279  occupancy_vs_ieta_HO->Fill(double(ieta));
280  }
281  //HF
282  if (sub == 4){
283  meSimHitsEnergyHF->Fill(en);
284  if (depth == 1){
285  emean_vs_ieta_HF1->Fill(double(ieta), en);
286  occupancy_vs_ieta_HF1->Fill(double(ieta));
287  }
288  if (depth == 2){
289  emean_vs_ieta_HF2->Fill(double(ieta), en);
290  occupancy_vs_ieta_HF2->Fill(double(ieta));
291  }
292  }
293  }
294 
295  //Ecal EB SimHits
297  ev.getByToken(tok_ecalEB_,ecalEBHits);
298  const PCaloHitContainer * SimHitResultEB = ecalEBHits.product () ;
299 
300  double EcalCone = 0;
301 
302  for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResultEB->begin () ; SimHits != SimHitResultEB->end(); ++SimHits) {
303 
304  EBDetId EBid = EBDetId(SimHits->id());
305 
306  const CaloCellGeometry* cellGeometry = geometry->getSubdetectorGeometry (EBid)->getGeometry (EBid) ;
307  double etaS = cellGeometry->getPosition().eta () ;
308  double phiS = cellGeometry->getPosition().phi () ;
309  double en = SimHits->energy();
310 
311  double r = dR(eta_MC, phi_MC, etaS, phiS);
312 
313  if (r < partR) EcalCone += en;
314  }
315 
316  //Ecal EE SimHits
318  ev.getByToken(tok_ecalEE_,ecalEEHits);
319  const PCaloHitContainer * SimHitResultEE = ecalEEHits.product () ;
320 
321  for (std::vector<PCaloHit>::const_iterator SimHits = SimHitResultEE->begin () ; SimHits != SimHitResultEE->end(); ++SimHits) {
322 
323  EEDetId EEid = EEDetId(SimHits->id());
324 
325  const CaloCellGeometry* cellGeometry = geometry->getSubdetectorGeometry (EEid)->getGeometry (EEid) ;
326  double etaS = cellGeometry->getPosition().eta () ;
327  double phiS = cellGeometry->getPosition().phi () ;
328  double en = SimHits->energy();
329 
330  double r = dR(eta_MC, phi_MC, etaS, phiS);
331 
332  if (r < partR) EcalCone += en;
333  }
334 
335  if (ietaMax != 0){ //If ietaMax == 0, there were no good HCAL SimHits
336  if (ietaMax > 0) ietaMax--; //Account for lack of ieta = 0
337 
338  meEnConeEtaProfile ->Fill(double(ietaMax), HcalCone);
339  meEnConeEtaProfile_E ->Fill(double(ietaMax), EcalCone);
340  meEnConeEtaProfile_EH ->Fill(double(ietaMax), HcalCone+EcalCone);
341  }
342 
343  nevtot++;
344 }
345 
346 
347 double HcalSimHitsValidation::dR(double eta1, double phi1, double eta2, double phi2) {
348  double PI = 3.1415926535898;
349  double deltaphi= phi1 - phi2;
350  if( phi2 > phi1 ) { deltaphi= phi2 - phi1;}
351  if(deltaphi > PI) { deltaphi = 2.*PI - deltaphi;}
352  double deltaeta = eta2 - eta1;
353  double tmp = sqrt(deltaeta* deltaeta + deltaphi*deltaphi);
354  return tmp;
355 }
356 
357 double HcalSimHitsValidation::phi12(double phi1, double en1, double phi2, double en2) {
358  // weighted mean value of phi1 and phi2
359 
360  double tmp;
361  double PI = 3.1415926535898;
362  double a1 = phi1; double a2 = phi2;
363 
364  if( a1 > 0.5*PI && a2 < 0.) a2 += 2*PI;
365  if( a2 > 0.5*PI && a1 < 0.) a1 += 2*PI;
366  tmp = (a1 * en1 + a2 * en2)/(en1 + en2);
367  if(tmp > PI) tmp -= 2.*PI;
368 
369  return tmp;
370 
371 }
372 
373 double HcalSimHitsValidation::dPhiWsign(double phi1, double phi2) {
374  // clockwise phi2 w.r.t phi1 means "+" phi distance
375  // anti-clockwise phi2 w.r.t phi1 means "-" phi distance
376 
377  double PI = 3.1415926535898;
378  double a1 = phi1; double a2 = phi2;
379  double tmp = a2 - a1;
380  if( a1*a2 < 0.) {
381  if(a1 > 0.5 * PI) tmp += 2.*PI;
382  if(a2 > 0.5 * PI) tmp -= 2.*PI;
383  }
384  return tmp;
385 
386 }
387 
388 
389 
391 
MonitorElement * meEnConeEtaProfile
T getUntrackedParameter(std::string const &, T const &) const
MonitorElement * meSimHitsEnergyHE
int i
Definition: DBlmapReader.cc:9
void setBinContent(int binx, double content)
set content of bin (1-D)
std::vector< PCaloHit > PCaloHitContainer
double dR(double eta1, double phi1, double eta2, double phi2)
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:873
#define PI
double phi12(double phi1, double en1, double phi2, double en2)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:434
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
MonitorElement * meSimHitsEnergyHO
double npart
Definition: HydjetWrapper.h:44
MonitorElement * emean_vs_ieta_HE1
edm::EDGetTokenT< edm::PCaloHitContainer > tok_ecalEB_
void Fill(long long x)
MonitorElement * meEnConeEtaProfile_E
MonitorElement * occupancy_vs_ieta_HB2
MonitorElement * occupancy_vs_ieta_HE2
edm::EDGetTokenT< edm::PCaloHitContainer > tok_ecalEE_
MonitorElement * meSimHitsEnergyHF
MonitorElement * emean_vs_ieta_HO
T sqrt(T t)
Definition: SSEVec.h:48
MonitorElement * occupancy_vs_ieta_HE3
void save(const std::string &filename, const std::string &path="", const std::string &pattern="", const std::string &rewrite="", const uint32_t run=0, SaveReferenceTag ref=SaveWithReference, int minStatus=dqm::qstatus::STATUS_OK, const std::string &fileupdate="RECREATE")
Definition: DQMStore.cc:2297
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:1187
edm::EDGetTokenT< edm::PCaloHitContainer > tok_hcal_
MonitorElement * get(const std::string &path) const
get ME from full pathname (e.g. &quot;my/long/dir/my_histo&quot;)
Definition: DQMStore.cc:1624
double dPhiWsign(double phi1, double phi2)
MonitorElement * emean_vs_ieta_HE2
tuple conf
Definition: dbtoconf.py:185
HcalSimHitsValidation(edm::ParameterSet const &conf)
MonitorElement * emean_vs_ieta_HE3
MonitorElement * emean_vs_ieta_HF1
MonitorElement * occupancy_vs_ieta_HB1
MonitorElement * emean_vs_ieta_HF2
const T & get() const
Definition: EventSetup.h:55
MonitorElement * meEnConeEtaProfile_EH
MonitorElement * emean_vs_ieta_HB1
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
double getBinContent(int binx) const
get content of bin (1-D)
MonitorElement * occupancy_vs_ieta_HF2
MonitorElement * occupancy_vs_ieta_HO
int getNbinsX(void) const
get # of bins in X-axis
edm::ESHandle< CaloGeometry > geometry
tuple cout
Definition: gather_cfg.py:121
MonitorElement * meSimHitsEnergyHB
edm::EDGetTokenT< edm::HepMCProduct > tok_evt_
virtual void analyze(edm::Event const &ev, edm::EventSetup const &c)
MonitorElement * emean_vs_ieta_HB2
MonitorElement * occupancy_vs_ieta_HF1
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:585
MonitorElement * occupancy_vs_ieta_HE1