CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HcalRecHitsMaker.cc
Go to the documentation of this file.
3 
15 #include "CLHEP/GenericFunctions/Erf.hh"
27 
28 #include "TFile.h"
29 #include "TGraph.h"
30 #include "TROOT.h"
31 #include <fstream>
32 #include <iostream>
33 
34 class RandomEngine;
35 
36 bool HcalRecHitsMaker::initialized_ = false;
41 std::vector<float> HcalRecHitsMaker::peds_;
42 std::vector<int> HcalRecHitsMaker::fctoadc_;
43 std::vector<float> HcalRecHitsMaker::sat_;
44 std::vector<float> HcalRecHitsMaker::gains_;
45 std::vector<float> HcalRecHitsMaker::noisesigma_;
46 std::vector<float> HcalRecHitsMaker::TPGFactor_;
47 std::vector<float> HcalRecHitsMaker::miscalib_;
48 std::vector<HcalDetId> HcalRecHitsMaker::theDetIds_;
49 std::vector<int> HcalRecHitsMaker::hbhi_;
50 std::vector<int> HcalRecHitsMaker::hehi_;
51 std::vector<int> HcalRecHitsMaker::hohi_;
52 std::vector<int> HcalRecHitsMaker::hfhi_;
53 unsigned HcalRecHitsMaker::maxIndex_ = 0 ;
54 
56  const RandomEngine * myrandom)
57  :
58  det_(det),
59  doDigis_(false),
60  random_(myrandom)
61 {
62  edm::ParameterSet RecHitsParameters=p.getParameter<edm::ParameterSet>("HCAL");
63  noise_ = RecHitsParameters.getParameter<std::vector<double> >("Noise");
64  corrfac_ = RecHitsParameters.getParameter<std::vector<double> >("NoiseCorrectionFactor");
65  threshold_ = RecHitsParameters.getParameter<std::vector<double> >("Threshold");
66  doSaturation_ = RecHitsParameters.getParameter<bool>("EnableSaturation");
67 
68  refactor_ = RecHitsParameters.getParameter<double> ("Refactor");
69  refactor_mean_ = RecHitsParameters.getParameter<double> ("Refactor_mean");
70  hcalfileinpath_= RecHitsParameters.getParameter<std::string> ("fileNameHcal");
71  inputCol_=RecHitsParameters.getParameter<edm::InputTag>("MixedSimHits");
73 
74  if(det_==4)
75  {
76  hbhi_.reserve(2600);
77  hehi_.reserve(2600);
78  }
79  else if (det_==5)
80  hohi_.reserve(2200);
81  else if (det_==6)
82  hfhi_.reserve(1800);
83 
84  nnoise_=noise_.size();
85  if(threshold_.size()!=noise_.size())
86  {
87  edm::LogWarning("CaloRecHitsProducer") << " WARNING : HCAL Noise simulation, the number of parameters should be the same for the noise and the thresholds. Disabling the noise simulation " << std::endl;
88  noise_.resize(nnoise_,0);
89  }
90 
91  // edm::ParameterSet hcalparam = p2.getParameter<edm::ParameterSet>("hcalSimParam");
92  // myHcalSimParameterMap_ = new HcalSimParameterMap(hcalparam);
93 
94  // Computes the fraction of HCAL above the threshold
95  Genfun::Erf myErf;
96  noiseFromDb_.resize(nnoise_,false);
97  hcalHotFraction_.resize(nnoise_,0.);
99  if(noise_.size()>0)
100  {
101  for(unsigned inoise=0;inoise<nnoise_;++inoise)
102  {
103  if(noise_[inoise]==0) {
104  hcalHotFraction_[inoise]=0.;
105  continue;
106  } else if(noise_[inoise]==-1) {
107  noiseFromDb_[inoise]=true;
108  continue;
109  } else {
110  hcalHotFraction_[inoise] = 0.5-0.5*myErf(threshold_[inoise]/noise_[inoise]/sqrt(2.));
111  myGaussianTailGenerators_[inoise]=new GaussianTail(random_,noise_[inoise],threshold_[inoise]);
112  }
113  }
114  }
115 }
116 
118 {
119  clean();
120  if(myGaussianTailGenerators_.size())
121  {
122  for(unsigned igt=0; igt<myGaussianTailGenerators_.size();++igt)
123  delete myGaussianTailGenerators_[igt];
124  }
126  theDetIds_.clear();
127  hbhi_.clear();
128  hehi_.clear();
129  hohi_.clear();
130  hfhi_.clear();
131 
132 }
133 
134 void HcalRecHitsMaker::init(const edm::EventSetup &es,bool doDigis,bool doMiscalib)
135 {
136  doDigis_=doDigis;
137  doMiscalib_=doMiscalib;
138 // needed for the noise simulation
139  edm::ESHandle<HcalDbService> conditions;
140  es.get<HcalDbRecord>().get(conditions);
141  const HcalDbService * theDbService=conditions.product();
142 
143  // get the correction factors
145  es.get<HcalRespCorrsRcd>().get(rchandle);
146  myRespCorr= rchandle.product();
147 
148  //general initialization
149  if(!initialized_)
150  {
151  theDetIds_.resize(9201);
152  unsigned ncells=createVectorsOfCells(es);
153  edm::LogInfo("CaloRecHitsProducer") << "Total number of cells in HCAL " << ncells << std::endl;
154  hcalRecHits_.resize(maxIndex_+1,0.);
155  edm::LogInfo("CaloRecHitsProducer") << "Largest HCAL hashedindex" << maxIndex_ << std::endl;
156 
157  peds_.resize(9201);
158  gains_.resize(9201);
159  sat_.resize(9201);
160  noisesigma_.resize(9201);
161 
162  miscalib_.resize(maxIndex_+1,1.);
163  // Read from file ( a la HcalRecHitsRecalib.cc)
164  // here read them from xml (particular to HCAL)
165  CaloMiscalibMapHcal mapHcal;
166 
168  es.get<IdealGeometryRecord>().get( topo );
169  mapHcal.prefillMap(*topo);
170 
171  edm::FileInPath hcalfiletmp("CalibCalorimetry/CaloMiscalibTools/data/"+hcalfileinpath_);
172  std::string hcalfile=hcalfiletmp.fullPath();
173  if(doMiscalib_ && !hcalfile.empty())
174  {
175  MiscalibReaderFromXMLHcal hcalreader_(mapHcal);
176 
177  hcalreader_.parseXMLMiscalibFile(hcalfile);
178  // mapHcal.print();
179  std::map<uint32_t,float>::const_iterator it=mapHcal.get().begin();
180  std::map<uint32_t,float>::const_iterator itend=mapHcal.get().end();
181  for(;it!=itend;++it)
182  {
183  HcalDetId myDetId(it->first);
184  float icalconst=it->second;
185  miscalib_[topo->detId2denseId(myDetId)]=refactor_mean_+(icalconst-refactor_mean_)*refactor_;
186  }
187  }
188 
189  // Open the histogram for the fC to ADC conversion
190  gROOT->cd();
191  edm::FileInPath myDataFile("FastSimulation/CaloRecHitsProducer/data/adcvsfc.root");
192  TFile * myFile = new TFile(myDataFile.fullPath().c_str(),"READ");
193  TGraph * myGraf = (TGraph*)myFile->Get("adcvsfc");
194  unsigned size=myGraf->GetN();
195  fctoadc_.resize(10000);
196  unsigned p_index=0;
197  fctoadc_[0]=0;
198  int prev_nadc=0;
199  int nadc=0;
200  for(unsigned ibin=0;ibin<size;++ibin)
201  {
202  double x,y;
203  myGraf->GetPoint(ibin,x,y);
204  int index=(int)floor(x);
205  if(index<0||index>=10000) continue;
206  prev_nadc=nadc;
207  nadc=(int)y;
208  // Now fills the vector
209  for(unsigned ivec=p_index;ivec<(unsigned)index;++ivec)
210  {
211  fctoadc_[ivec] = prev_nadc;
212  }
213  p_index = index;
214  }
215  myFile->Close();
216  gROOT->cd();
217  edm::FileInPath myTPGFilePath("CalibCalorimetry/HcalTPGAlgos/data/RecHit-TPG-calib.dat");
218  TPGFactor_.resize(87,1.2);
219  std::ifstream myTPGFile(myTPGFilePath.fullPath().c_str(),ifstream::in);
220  if(myTPGFile)
221  {
222  float gain;
223  myTPGFile >> gain;
224  for(unsigned i=0;i<86;++i)
225  {
226  myTPGFile >> TPGFactor_[i] ;
227  //std::cout << TPGFactor_[i] << std::endl;
228  }
229  }
230  else
231  {
232  std::cout << " Unable to open CalibCalorimetry/HcalTPGAlgos/data/RecHit-TPG-calib.dat" << std::endl;
233  std::cout << " Using a constant 1.2 factor " << std::endl;
234  }
235 
236  initialized_=true;
237  }
238 
239  if(!initializedHB_ && det_==4){
240  nhbcells_ = hbhi_.size();
241  //HB
242  for(unsigned ic=0;ic<nhbcells_;++ic)
243  {
244  float gain = theDbService->getGain(theDetIds_[hbhi_[ic]])->getValue(0);
245  float mgain=0.;
246  for(unsigned ig=0;ig<4;++ig)
247  mgain+=theDbService->getGain(theDetIds_[hbhi_[ic]])->getValue(ig);
248  if(noiseFromDb_[0])
249  noisesigma_[hbhi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hbhi_[ic]])*mgain*0.25;
250  // std::cout << " NOISEHB " << theDetIds_[hbhi_[ic]].ieta() << " " << noisesigma_[hbhi_[ic]] << " "<< std::endl;
251  // 10000 (ADC scale) / 4. (to compute the mean) / 0.92 ADC/fC
252  // *1.25 (only ~80% in 1ts Digi, while saturation applied to 4ts RecHit)
253  int ieta=theDetIds_[hbhi_[ic]].ieta();
254  float tpgfactor=TPGFactor_[(ieta>0)?ieta+43:-ieta];
255  mgain*=2500./0.92*tpgfactor ;// 10000 (ADC scale) / 4. (to compute the mean)
256  sat_[hbhi_[ic]]=(doSaturation_)?mgain:99999.;
257 
258  peds_[hbhi_[ic]]=theDbService->getPedestal(theDetIds_[hbhi_[ic]])->getValue(0);
259 
260  gain*=tpgfactor;
261  gains_[hbhi_[ic]]=gain;
262  }
263 
264  initializedHB_=true;
265  }
266 
267  if(!initializedHE_ && det_==4){
268  nhecells_ = hehi_.size();
269  //HE
270  for(unsigned ic=0;ic<nhecells_;++ic)
271  {
272  float gain= theDbService->getGain(theDetIds_[hehi_[ic]])->getValue(0);
273  int ieta=theDetIds_[hehi_[ic]].ieta();
274  float mgain=0.;
275  for(unsigned ig=0;ig<4;++ig)
276  mgain+=theDbService->getGain(theDetIds_[hehi_[ic]])->getValue(ig);
277  if(noiseFromDb_[1])
278  noisesigma_[hehi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hehi_[ic]])*mgain*0.25;
279  // std::cout << " NOISEHE " << theDetIds_[hehi_[ic]].ieta() << " " << noisesigma_[hehi_[ic]] << " "<< std::endl;
280  float tpgfactor=TPGFactor_[(ieta>0)?ieta+44:-ieta+1];
281  mgain*=2500./0.92*tpgfactor;
282  sat_[hehi_[ic]]=(doSaturation_)?mgain:99999.;
283 
284  gain*=tpgfactor;
285  peds_[hehi_[ic]]=theDbService->getPedestal(theDetIds_[hehi_[ic]])->getValue(0);
286  gains_[hehi_[ic]]=gain;
287  }
288 
289  initializedHE_=true;
290  }
291 
292  if(!initializedHO_ && det_==5){
293  nhocells_ = hohi_.size();
294  //HO
295  for(unsigned ic=0;ic<nhocells_;++ic)
296  {
297  float ped=theDbService->getPedestal(theDetIds_[hohi_[ic]])->getValue(0);
298  float gain=theDbService->getGain(theDetIds_[hohi_[ic]])->getValue(0);
299  float mgain=0.;
300  for(unsigned ig=0;ig<4;++ig)
301  mgain+=theDbService->getGain(theDetIds_[hohi_[ic]])->getValue(ig);
302  if(noiseFromDb_[0])
303  noisesigma_[hohi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hohi_[ic]])*mgain*0.25;
304  // std::cout << " NOISEHO " << theDetIds_[hohi_[ic]].ieta() << " " << noisesigma_[hohi_[ic]] << " "<< std::endl;
305  int ieta=HcalDetId(hohi_[ic]).ieta();
306  float tpgfactor=TPGFactor_[(ieta>0)?ieta+43:-ieta];
307  mgain*=2500./0.92*tpgfactor;
308  sat_[hohi_[ic]]=(doSaturation_)?mgain:99999.;
309 
310  gain*=tpgfactor;
311  peds_[hohi_[ic]]=ped;
312  gains_[hohi_[ic]]=gain;
313  }
314 
315  initializedHO_=true;
316  }
317 
318  if(!initializedHF_ && det_==6){
319  nhfcells_ = hfhi_.size();
320  //HF
321  for(unsigned ic=0;ic<nhfcells_;++ic)
322  {
323  float ped=theDbService->getPedestal(theDetIds_[hfhi_[ic]])->getValue(0);
324  float gain=theDbService->getGain(theDetIds_[hfhi_[ic]])->getValue(0);
325  float mgain=0.;
326  for(unsigned ig=0;ig<4;++ig)
327  mgain+=theDbService->getGain(theDetIds_[hfhi_[ic]])->getValue(ig);
328  // additional 1/2 factor for the HF (Salavat)
329  if(noiseFromDb_[0])
330  {
331  // computation from when the noise was taken
332  noisesigma_[hfhi_[ic]]=noiseInfCfromDB(theDbService,theDetIds_[hfhi_[ic]])*mgain*0.25;
333  }
334  // std::cout << " NOISEHF " << theDetIds_[hfhi_[ic]].ieta() << " " << noisesigma_[hfhi_[ic]] << " "<< std::endl;
335 
336  mgain*=2500./0.36;
337  sat_[hfhi_[ic]]=(doSaturation_)?mgain:99999.;
338  int ieta=theDetIds_[hfhi_[ic]].ieta();
339  gain*=TPGFactor_[(ieta>0)?ieta+45:-ieta+2];
340  peds_[hfhi_[ic]]=ped;
341  gains_[hfhi_[ic]]=gain;
342  }
343 
344  initializedHF_=true;
345  }
346 
347  // clear the vector we don't need. It is a bit stupid
348  hcalRecHits_.resize(maxIndex_+1,0.);
349 
350 }
351 
352 
353 // Get the PCaloHits from the event. They have to be stored in a map, because when
354 // the pile-up is added thanks to the Mixing Module, the same cell can be present several times
356 {
357  clean();
358 
360  iEvent.getByLabel(inputCol_,cf);
361  std::auto_ptr<MixCollection<PCaloHit> > colcalo(new MixCollection<PCaloHit>(cf.product(),std::pair<int,int>(0,0) ));
362 
363  MixCollection<PCaloHit>::iterator it=colcalo->begin();;
364  MixCollection<PCaloHit>::iterator itend=colcalo->end();
365  unsigned counter=0;
366  for(;it!=itend;++it)
367  {
368  HcalDetId detid(it->id());
369  int hashedindex=topo.detId2denseId(detid);
370 
371  // apply ToF correction
372  int time_slice=0; // temporary
373  double fTOF=1.;
374  if (detid.subdet()==HcalForward) fTOF = (time_slice==0) ? 1. : 0.;
375  else fTOF = fractionOOT(time_slice);
376 
377  switch(detid.subdet())
378  {
379  case HcalBarrel:
380  {
381  if(det_==4)
382  Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[0],corrfac_[0]);
383  }
384  break;
385  case HcalEndcap:
386  {
387  if(det_==4)
388  Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[1],corrfac_[1]);
389  }
390  break;
391  case HcalOuter:
392  {
393  if(det_==5)
394  Fill(hashedindex,fTOF*(it->energy()),firedCells_,noise_[0],corrfac_[0]);
395  }
396  break;
397  case HcalForward:
398  {
399  if(det_==6 && time_slice==0) // skip the HF hit if out-of-time
400  Fill(hashedindex,it->energy(),firedCells_,noise_[0],corrfac_[0]);
401  }
402  break;
403  default:
404  edm::LogWarning("CaloRecHitsProducer") << "RecHit not registered\n";
405  ;
406  }
407  ++counter;
408  }
409 }
410 
411 // Fills the collections.
413 {
414  loadPCaloHits(iEvent,topo);
415  noisify();
416  hbheHits.reserve(firedCells_.size());
417  if(doDigis_)
418  {
419  hbheDigis.reserve(firedCells_.size());
420  }
421  static HcalQIESample zeroSample(0,0,0,0);
422  unsigned nhits=firedCells_.size();
423  // HB and HE
424 
425  for(unsigned ihit=0;ihit<nhits;++ihit)
426  {
427  unsigned cellhashedindex=firedCells_[ihit];
428  const HcalDetId& detid = theDetIds_[cellhashedindex];
429  unsigned subdet=(detid.subdet()==HcalBarrel) ? 0: 1;
430 
431  float energy=hcalRecHits_[cellhashedindex];
432  // Check if it is above the threshold
433  if(energy<threshold_[subdet]) continue;
434  // apply RespCorr only to the RecHit
435  energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
436  // poor man saturation
437  if(energy>sat_[cellhashedindex])
438  {
439  // std::cout << " Saturation " << energy << " " << sat_[cellhashedindex] << " HB " << std::endl;
440  energy=sat_[cellhashedindex];
441  }
442 
443 
444  hbheHits.push_back(HBHERecHit(detid,energy,0.));
445  if(doDigis_)
446  {
447  HBHEDataFrame myDataFrame(detid);
448  myDataFrame.setSize(2);
449  double nfc=hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
450  int nadc=fCtoAdc(nfc);
451  HcalQIESample qie(nadc, 0, 0, 0) ;
452  myDataFrame.setSample(0,qie);
453  myDataFrame.setSample(1,zeroSample);
454  hbheDigis.push_back(myDataFrame);
455  }
456  }
457 }
458 
459 
460 // Fills the collections.
462 {
463  loadPCaloHits(iEvent,topo);
464  noisify();
465  hoHits.reserve(firedCells_.size());
466  if(doDigis_)
467  {
468  hoDigis.reserve(firedCells_.size());
469  }
470  static HcalQIESample zeroSample(0,0,0,0);
471 
472  // HO
473  unsigned nhits=firedCells_.size();
474  for(unsigned ihit=0;ihit<nhits;++ihit)
475  {
476 
477  unsigned cellhashedindex=firedCells_[ihit];
478  // Check if it is above the threshold
479  if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue;
480 
481  const HcalDetId& detid=theDetIds_[cellhashedindex];
482  int ieta = detid.ieta();
483 
484  // Force only Ring#0 to remain
485  if(ieta > 4 || ieta < -4 ) continue;
486 
487  float energy=hcalRecHits_[cellhashedindex];
488  // apply RespCorr
489  energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
490 
491  // poor man saturation
492  if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
493 
494  hoHits.push_back(HORecHit(detid,energy,0));
495  }
496 }
497 
498 // Fills the collections.
500 {
501  loadPCaloHits(iEvent,topo);
502  noisify();
503  hfHits.reserve(firedCells_.size());
504  if(doDigis_)
505  {
506  hfDigis.reserve(firedCells_.size());
507  }
508  static HcalQIESample zeroSample(0,0,0,0);
509 
510  unsigned nhits=firedCells_.size();
511  for(unsigned ihit=0;ihit<nhits;++ihit)
512  {
513  unsigned cellhashedindex=firedCells_[ihit];
514  // Check if it is above the threshold
515  if(hcalRecHits_[cellhashedindex]<threshold_[0]) continue;
516 
517  float energy=hcalRecHits_[cellhashedindex];
518 
519  // apply RespCorr
520  energy *= myRespCorr->getValues(theDetIds_[cellhashedindex])->getValue();
521 
522  // poor man saturation
523  if(energy>sat_[cellhashedindex]) energy=sat_[cellhashedindex];
524 
525  const HcalDetId & detid=theDetIds_[cellhashedindex];
526  hfHits.push_back(HFRecHit(detid,energy,0.));
527  if(doDigis_)
528  {
529  HFDataFrame myDataFrame(detid);
530  myDataFrame.setSize(1);
531 
532  double nfc= hcalRecHits_[cellhashedindex]/gains_[cellhashedindex]+peds_[cellhashedindex];
533  int nadc=fCtoAdc(nfc/2.6);
534  HcalQIESample qie(nadc, 0, 0, 0) ;
535  myDataFrame.setSample(0,qie);
536  hfDigis.push_back(myDataFrame);
537  }
538  }
539 }
540 
541 
542 // For a fast injection of the noise: the list of cell ids is stored
544 {
546  es.get<CaloGeometryRecord>().get(pG);
548  es.get<IdealGeometryRecord>().get( topo );
549 
554 
556 }
557 
558 // list of the cellids for a given subdetector
559 unsigned HcalRecHitsMaker::createVectorOfSubdetectorCells(const CaloGeometry& cg, const HcalTopology& topo, int subdetn,std::vector<int>& cellsvec )
560 {
561 
563  const std::vector<DetId>& ids=geom->getValidDetIds(DetId::Hcal,subdetn);
564 
565  for (std::vector<DetId>::const_iterator i=ids.begin(); i!=ids.end(); i++)
566  {
567  HcalDetId myDetId(*i);
568  // std::cout << myDetId << myHcalSimParameterMap_->simParameters(myDetId).simHitToPhotoelectrons() << std::endl;;
569  // std::cout << " hi " << hi << " " theDetIds_.size() << std::endl;
570  unsigned int hi=topo.detId2denseId(myDetId);
571  theDetIds_[hi]=myDetId;
572  // std::cout << myDetId << " " << hi << std::endl;
573  cellsvec.push_back(hi);
574 
575  if(hi>maxIndex_)
576  maxIndex_=hi;
577  }
578  return cellsvec.size();
579 }
580 
581 // Takes a hit (from a PSimHit) and fills a map
582 void HcalRecHitsMaker::Fill(int id, float energy, std::vector<int>& theHits,float noise,float correctionfactor)
583 {
584  if(doMiscalib_)
585  energy*=miscalib_[id];
586 
587  if(noise==-1)
588  noise=noisesigma_[id]*correctionfactor;
589 
590  // Check if the RecHit exists
591  if(hcalRecHits_[id]>0.)
592  hcalRecHits_[id]+=energy;
593  else
594  {
595  // the noise is injected only the first time
596  hcalRecHits_[id]=energy + random_->gaussShoot(0.,noise);
597  theHits.push_back(id);
598  }
599 }
600 
602 {
603  unsigned total=0;
604  switch(det_)
605  {
606  case 4:
607  {
608  // do the HB
609  if(noise_[0] != 0.) {
611  }
612  // do the HE
613  if(noise_[1] != 0.) {
615  }
616  }
617  break;
618  case 5:
619  {
620  // do the HO
621  if(noise_[0] != 0.) {
623  }
624  }
625  break;
626  case 6:
627  {
628  // do the HF
629  if(noise_[0] != 0.) {
631  }
632  }
633  break;
634  default:
635  break;
636  }
637  edm::LogInfo("CaloRecHitsProducer") << "CaloRecHitsProducer : added noise in "<< total << " HCAL cells " << std::endl;
638 }
639 
640 unsigned HcalRecHitsMaker::noisifySubdet(std::vector<float>& theMap, std::vector<int>& theHits, const std::vector<int>& thecells, unsigned ncells, double hcalHotFraction,const GaussianTail *myGT,double noise_,double threshold,double correctionfactor)
641 {
642  // If the fraction of "hot " is small, use an optimized method to inject noise only in noisy cells. The 30% has not been tuned
643  if(noise_!=-1 && hcalHotFraction==0.) return 0;
644  if(hcalHotFraction<0.3 && noise_!=-1)
645  {
646  double mean = (double)(ncells-theHits.size())*hcalHotFraction;
647  unsigned nhcal = random_->poissonShoot(mean);
648 
649  unsigned ncell=0;
650  unsigned cellindex=0;
651  uint32_t cellhashedindex=0;
652 
653  while(ncell < nhcal)
654  {
655  cellindex = (unsigned)(random_->flatShoot()*ncells);
656  cellhashedindex = thecells[cellindex];
657 
658  if(hcalRecHits_[cellhashedindex]==0.) // new cell
659  {
660  hcalRecHits_[cellhashedindex]=myGT->shoot();
661  theHits.push_back(cellhashedindex);
662  ++ncell;
663  }
664  }
665  return ncell;
666  }
667  else // otherwise, inject noise everywhere
668  {
669  uint32_t cellhashedindex=0;
670  unsigned nhcal=thecells.size();
671 
672  double sigma = noise_;
673  for(unsigned ncell=0;ncell<nhcal;++ncell)
674  {
675  cellhashedindex = thecells[ncell];
676  if(hcalRecHits_[cellhashedindex]==0.) // new cell
677  {
678  //only use noisesigma if it has been filled from the DB
679  //otherwise, leave sigma alone (default is user input value)
680  if(noise_==-1) sigma=noisesigma_[cellhashedindex]*correctionfactor;
681 
682  double noise =random_->gaussShoot(0.,sigma);
683  if(noise>threshold)
684  {
685  hcalRecHits_[cellhashedindex]=noise;
686  theHits.push_back(cellhashedindex);
687  }
688  }
689  }
690  return nhcal;
691  }
692  return 0;
693 }
694 
696 {
698 }
699 
700 void HcalRecHitsMaker::cleanSubDet(std::vector<float>& hits,std::vector<int>& cells)
701 {
702  unsigned size=cells.size();
703  // Reset the energies
704  for(unsigned ic=0;ic<size;++ic)
705  {
706  hits[cells[ic]] = 0.;
707  }
708  // Clear the list of fired cells
709  cells.clear();
710 }
711 
712 // fC to ADC conversion
713 int HcalRecHitsMaker::fCtoAdc(double fc) const
714 {
715  if(fc<0.) return 0;
716  if(fc>9985.) return 127;
717  return fctoadc_[(unsigned)floor(fc)];
718 }
719 
720 double HcalRecHitsMaker::noiseInfCfromDB(const HcalDbService * conditions,const HcalDetId & detId)
721 {
722  // method from Salavat
723  // fetch pedestal widths (noise sigmas for all 4 CapID)
724  const HcalPedestalWidth* pedWidth =
725  conditions-> getPedestalWidth(detId);
726  double ssqq_1 = pedWidth->getSigma(0,0);
727  double ssqq_2 = pedWidth->getSigma(1,1);
728  double ssqq_3 = pedWidth->getSigma(2,2);
729  double ssqq_4 = pedWidth->getSigma(3,3);
730 
731  // correction factors (hb,he,ho,hf)
732  // static float corrfac[4]={1.39,1.32,1.17,3.76};
733  // static float corrfac[4]={0.93,0.88,1.17,2.51}; // divided HB, HE, HF by 1.5 (to take into account the halving of the number of time slices), HO did not change
734 
735  int sub = detId.subdet();
736 
737  // HO: only Ring#0 matters
738  int ieta = detId.ieta();
739  if (sub == 3 && abs (ieta) > 4) return 0.;
740 
741  // effective RecHits (for this particular detId) noise calculation :
742 
743  double sig_sq_mean = 0.25 * ( ssqq_1 + ssqq_2 + ssqq_3 + ssqq_4);
744 
745  // f - external parameter (currently set to 0.5 in the FullSim) !!!
746  double f=0.5;
747 
748  double term = sqrt (1. + sqrt(1. - 2.*f*f));
749  double alpha = sqrt (0.5) * term;
750  double beta = sqrt (0.5) * f / term;
751 
752  // double RMS1 = sqrt(sig_sq_mean) * sqrt (2.*beta*beta + alpha*alpha) ;
753  double RMS4 = sqrt(sig_sq_mean) * sqrt (2.*beta*beta + 2.*(alpha-beta)*(alpha-beta) + 2.*(alpha-2.*beta)*(alpha-2.*beta)) ;
754 
755  double noise_rms_fC;
756 
757  // if(sub == 4) noise_rms_fC = RMS1;
758  // else noise_rms_fC = RMS4;
759  noise_rms_fC = RMS4;
760 
761 
762  // to convert from above fC to GeV - multiply by gain (GeV/fC)
763  // const HcalGain* gain = conditions->getGain(detId);
764  // double noise_rms_GeV = noise_rms_fC * gain->getValue(0); // Noise RMS (GeV) for detId channel
765  return noise_rms_fC;
766 }
767 
768 // fraction of energy collected as a function of ToF (for out-of-time particles; use case is out-of-time pileup)
769 double HcalRecHitsMaker::fractionOOT(int time_slice)// in units of 25 ns; 0 means in-time
770 {
771  double SF = 100./88.; // to normalize to in-time signal (88% is in the reco window)
772  if (time_slice==-4) {
773  return 0.02*SF;
774  } else if (time_slice==-3) {
775  return 0.06*SF;
776  } else if (time_slice==-2) {
777  return 0.19*SF;
778  } else if (time_slice==-1) {
779  return 0.24*SF;
780  } else if (time_slice==0) {
781  return 1.;
782  } else if (time_slice==1) {
783  return 0.70*SF;
784  } else {
785  return 0.;
786  }
787 
788 }
const double beta
static const char noise_[]
T getParameter(std::string const &) const
const CaloSubdetectorGeometry * getSubdetectorGeometry(const DetId &id) const
access the subdetector geometry for the given subdetector directly
Definition: CaloGeometry.cc:43
const std::map< uint32_t, float > & get()
int i
Definition: DBlmapReader.cc:9
static bool initializedHO_
static std::vector< float > peds_
float alpha
Definition: AMPTWrapper.h:95
HcalRecHitsMaker(edm::ParameterSet const &p, int, const RandomEngine *random)
virtual unsigned int detId2denseId(const DetId &id) const
return a linear packed id
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:32
static std::vector< float > sat_
std::vector< GaussianTail * > myGaussianTailGenerators_
void Fill(int id, float energy, std::vector< int > &myHits, float noise, float correctionfactor)
float getSigma(int fCapId1, int fCapId2) const
get correlation element for capId1/2 = 0..3
int fCtoAdc(double fc) const
std::vector< double > hcalHotFraction_
double noiseInfCfromDB(const HcalDbService *conditions, const HcalDetId &detId)
edm::InputTag inputCol_
std::vector< double > noise_
bool parseXMLMiscalibFile(std::string configFile)
std::vector< float > hcalRecHits_
static std::vector< float > TPGFactor_
void push_back(T const &t)
#define abs(x)
Definition: mlp_lapack.h:159
const Item * getValues(DetId fId, bool throwOnFail=true) const
const RandomEngine * random_
static std::vector< int > hfhi_
float getValue(int fCapId) const
get value for capId = 0..3
Definition: HcalGain.h:20
unsigned createVectorsOfCells(const edm::EventSetup &es)
void setSize(int size)
std::vector< double > corrfac_
double gaussShoot(double mean=0.0, double sigma=1.0) const
Definition: RandomEngine.h:37
virtual const std::vector< DetId > & getValidDetIds(DetId::Detector det=DetId::Detector(0), int subdet=0) const
Get a list of valid detector ids (for the given subdetector)
static std::vector< float > noisesigma_
double fractionOOT(int time_slice)
void loadHcalRecHits(edm::Event &iEvent, const HcalTopology &, HBHERecHitCollection &hbheHits, HBHEDigiCollection &hbheDigis)
int iEvent
Definition: GenABIO.cc:243
std::vector< double > threshold_
static std::vector< int > hohi_
static std::vector< int > fctoadc_
void setSample(int i, const HcalQIESample &sam)
Definition: HBHEDataFrame.h:50
T sqrt(T t)
Definition: SSEVec.h:48
unsigned int poissonShoot(double mean) const
Definition: RandomEngine.h:44
static bool initializedHE_
int ieta() const
get the cell ieta
Definition: HcalDetId.h:38
static bool initialized_
float getValue(int fCapId) const
get value for capId = 0..3
Definition: HcalPedestal.h:19
double f[11][100]
std::string hcalfileinpath_
static std::vector< HcalDetId > theDetIds_
void setSample(int i, const HcalQIESample &sam)
Definition: HFDataFrame.h:50
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
unsigned createVectorOfSubdetectorCells(const CaloGeometry &, const HcalTopology &, int subdetn, std::vector< int > &)
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:62
const HcalGain * getGain(const HcalGenericDetId &fId) const
T const * product() const
Definition: Handle.h:74
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngine.h:30
static unsigned maxIndex_
void setSize(int size)
Definition: HFDataFrame.cc:17
const HcalRespCorrs * myRespCorr
void prefillMap(const HcalTopology &topology)
unsigned noisifySubdet(std::vector< float > &theMap, std::vector< int > &theHits, const std::vector< int > &thecells, unsigned ncells, double hcalHotFraction_, const GaussianTail *, double sigma, double threshold, double correctionfactor)
double shoot() const
Definition: GaussianTail.cc:20
static std::vector< float > gains_
float getValue() const
Definition: HcalRespCorr.h:18
void reserve(size_type n)
static bool initializedHB_
void cleanSubDet(std::vector< float > &hits, std::vector< int > &cells)
std::vector< int > firedCells_
static std::vector< int > hbhi_
tuple cout
Definition: gather_cfg.py:121
void loadPCaloHits(const edm::Event &iEvent, const HcalTopology &)
std::string fullPath() const
Definition: FileInPath.cc:171
Definition: DDAxes.h:10
static std::vector< float > miscalib_
std::vector< bool > noiseFromDb_
const HcalPedestal * getPedestal(const HcalGenericDetId &fId) const
static bool initializedHF_
static std::vector< int > hehi_
void init(const edm::EventSetup &es, bool dodigis, bool domiscalib)
tuple size
Write out results.