CMS 3D CMS Logo

CalorimetryManager.cc
Go to the documentation of this file.
1 //updated by Reza Goldouzian
2 //Framework headers
5 
6 // Fast Simulation headers
18 //#include "FastSimulation/Utilities/interface/Histos.h"
28 // New headers for Muon Mip Simulation
31 // Muon Brem
33 
34 //Gflash Hadronic Model
42 
43 //FastHFShowerLibrary
45 
46 // STL headers
47 #include <vector>
48 #include <iostream>
49 
50 //CMSSW headers
53 //#include "DataFormats/EcalDetId/interface/EcalDetId.h"
55 
58 
59 //ROOT headers
60 #include "TROOT.h"
61 #include "TH1.h"
62 
63 using namespace edm;
64 
67 
68 std::vector<std::pair<int,float> > CalorimetryManager::myZero_ = std::vector<std::pair<int,float> >
69  (1,std::pair<int,float>(0,0.));
70 
72  myCalorimeter_(nullptr),
73  // myHistos(0),
74  initialized_(false)
75 {;}
76 
78  const edm::ParameterSet& fastCalo,
79  const edm::ParameterSet& fastMuECAL,
80  const edm::ParameterSet& fastMuHCAL,
81  const edm::ParameterSet& parGflash)
82  :
83  mySimEvent(aSimEvent),
86 {
87 
90 
91  //Gflash
92  theProfile = new GflashHadronShowerProfile(parGflash);
93  thePiKProfile = new GflashPiKShowerProfile(parGflash);
96 
97  // FastHFShowerLibrary
99 
100  readParameters(fastCalo);
101 
102  // myHistos = 0;
103 
104  // myHistos = Histos::instance();
105  // myHistos->book("h10",140,-3.5,3.5,100,-0.5,99.5);
106  // myHistos->book("h20",150,0,150.,100,-0.5,99.5);
107  // myHistos->book("h100",140,-3.5,3.5,100,0,0.1);
108  // myHistos->book("h110",140,-3.5,3.5,100,0,10.);
109  // myHistos->book("h120",200,-5.,5.,100,0,0.5);
110 
111  // myHistos->book("h200",300,0,3.,100,0.,35.);
112  // myHistos->book("h210",720,-M_PI,M_PI,100,0,35.);
113  // myHistos->book("h212",720,-M_PI,M_PI,100,0,35.);
114 
115  // myHistos->bookByNumber("h30",0,7,300,-3.,3.,100,0.,35.);
116  // myHistos->book("h310",75,-3.,3.,"");
117  // myHistos->book("h400",100,-10.,10.,100,0.,35.);
118  // myHistos->book("h410",720,-M_PI,M_PI);
119 
120  myCalorimeter_ =
121  new CaloGeometryHelper(fastCalo);
122  myHDResponse_ =
123  new HCALResponse(fastCalo.getParameter<edm::ParameterSet>("HCALResponse"));
124  myHSParameters_ =
125  new HSParameters(fastCalo.getParameter<edm::ParameterSet>("HSParameters"));
126 
127  // Material Effects for Muons in ECAL (only EnergyLoss implemented so far)
128 
129  if ( fastMuECAL.getParameter<bool>("PairProduction") ||
130  fastMuECAL.getParameter<bool>("Bremsstrahlung") ||
131  fastMuECAL.getParameter<bool>("MuonBremsstrahlung") ||
132  fastMuECAL.getParameter<bool>("EnergyLoss") ||
133  fastMuECAL.getParameter<bool>("MultipleScattering") )
134  theMuonEcalEffects = new MaterialEffects(fastMuECAL);
135 
136  // Material Effects for Muons in HCAL (only EnergyLoss implemented so far)
137 
138  if ( fastMuHCAL.getParameter<bool>("PairProduction") ||
139  fastMuHCAL.getParameter<bool>("Bremsstrahlung") ||
140  fastMuHCAL.getParameter<bool>("MuonBremsstrahlung") ||
141  fastMuHCAL.getParameter<bool>("EnergyLoss") ||
142  fastMuHCAL.getParameter<bool>("MultipleScattering") )
143  theMuonHcalEffects = new MaterialEffects(fastMuHCAL);
144 
145  if( fastCalo.exists("ECALResponseScaling") ) {
146  ecalCorrection = std::unique_ptr<KKCorrectionFactors>( new KKCorrectionFactors( fastCalo.getParameter<edm::ParameterSet>("ECALResponseScaling") ) );
147  }
148 
149 }
150 
152 {
153  EBMapping_.clear();
154  EEMapping_.clear();
155  HMapping_.clear();
156  ESMapping_.clear();
157  muonSimTracks.clear();
158 }
159 
161 {
162  if(myCalorimeter_) delete myCalorimeter_;
163  if(myHDResponse_) delete myHDResponse_;
164 
167 
168  if ( theProfile ) delete theProfile;
169 
171 }
172 
174 {
175 
176  if(!evtsToDebug_.empty())
177  {
178  std::vector<unsigned int>::const_iterator itcheck=find(evtsToDebug_.begin(),evtsToDebug_.end(),mySimEvent->id().event());
179  debug_=(itcheck!=evtsToDebug_.end());
180  if(debug_)
181  mySimEvent->print();
182  }
183  // Clear the content of the calorimeters
184  if(!initialized_)
185  {
186 
187  // Check if the preshower is really available
189  {
190  std::cout << " WARNING " << std::endl;
191  std::cout << " The preshower simulation has been turned on; but no preshower geometry is available " << std::endl;
192  std::cout << " Disabling the preshower simulation " << std::endl;
193  simulatePreshower_ = false;
194  }
195 
196  initialized_=true;
197  }
198  clean();
199 
200  LogInfo("FastCalorimetry") << "Reconstructing " << (int) mySimEvent->nTracks() << " tracks." << std::endl;
201  for( int fsimi=0; fsimi < (int) mySimEvent->nTracks() ; ++fsimi) {
202 
203  FSimTrack& myTrack = mySimEvent->track(fsimi);
204 
205  int pid = abs(myTrack.type());
206 
207  if (debug_) {
208  LogInfo("FastCalorimetry") << " ===> pid = " << pid << std::endl;
209  }
210 
211 
212  // Check that the particle hasn't decayed
213  if(myTrack.noEndVertex()) {
214  // Simulate energy smearing for photon and electrons
215  if ( pid == 11 || pid == 22 ) {
216 
217  if ( myTrack.onEcal() )
218  EMShowerSimulation(myTrack, random);
219  else if ( myTrack.onVFcal() ) {
220  if(useShowerLibrary) {
222  myHDResponse_->correctHF(myTrack.hcalEntrance().e(),abs(myTrack.type()));
224  }
225  else reconstructHCAL(myTrack, random);
226  }
227  } // electron or photon
228  else if (pid==13)
229  {
230  MuonMipSimulation(myTrack, random);
231  }
232  // Simulate energy smearing for hadrons (i.e., everything
233  // but muons... and SUSY particles that deserve a special
234  // treatment.
235  else if ( pid < 1000000 ) {
236  if ( myTrack.onHcal() || myTrack.onVFcal() ) {
237  if(optionHDSim_ == 0 ) reconstructHCAL(myTrack, random);
238  else HDShowerSimulation(myTrack, random);
239  }
240  } // pid < 1000000
241  } // myTrack.noEndVertex()
242  } // particle loop
243  // LogInfo("FastCalorimetry") << " Number of hits (barrel)" << EBMapping_.size() << std::endl;
244  // LogInfo("FastCalorimetry") << " Number of hits (Hcal)" << HMapping_.size() << std::endl;
245  // std::cout << " Nombre de hit (endcap)" << EEMapping_.size() << std::endl;
246 
247 } // reconstruct
248 
249 // Simulation of electromagnetic showers in PS, ECAL, HCAL
252  std::vector<const RawParticle*> thePart;
253  double X0depth;
254 
255  if (debug_) {
256  LogInfo("FastCalorimetry") << " EMShowerSimulation " <<myTrack << std::endl;
257  }
258 
259  // std::cout << " Simulating " << myTrack << std::endl;
260 
261  // The Particle at ECAL entrance
262  // std::cout << " Before ecalEntrance " << std::endl;
263  myPart = myTrack.ecalEntrance();
264 
265  // protection against infinite loop.
266  if ( myTrack.type() == 22 && myPart.e()<0.055) return;
267 
268 
269  // Barrel or Endcap ?
270  int onEcal = myTrack.onEcal();
271  int onHcal = myTrack.onHcal();
272  int onLayer1 = myTrack.onLayer1();
273  int onLayer2 = myTrack.onLayer2();
274 
275  // The entrance in ECAL
276  XYZPoint ecalentrance = myPart.vertex().Vect();
277 
278  // std::cout << " Ecal entrance " << ecalentrance << std::endl;
279 
280  // The preshower
281  PreshowerHitMaker * myPreshower = nullptr;
282  if(simulatePreshower_ && (onLayer1 || onLayer2))
283  {
284  XYZPoint layer1entrance,layer2entrance;
285  XYZVector dir1,dir2;
286  if(onLayer1)
287  {
288  layer1entrance = XYZPoint(myTrack.layer1Entrance().vertex().Vect());
289  dir1 = XYZVector(myTrack.layer1Entrance().Vect().Unit());
290  }
291  if(onLayer2)
292  {
293  layer2entrance = XYZPoint(myTrack.layer2Entrance().vertex().Vect());
294  dir2 = XYZVector(myTrack.layer2Entrance().Vect().Unit());
295  }
296  // std::cout << " Layer1entrance " << layer1entrance << std::endl;
297  // std::cout << " Layer2entrance " << layer2entrance << std::endl;
298  myPreshower = new PreshowerHitMaker(myCalorimeter_,
299  layer1entrance,
300  dir1,
301  layer2entrance,
302  dir2,
304  random);
305  myPreshower->setMipEnergy(mipValues_[0],mipValues_[1]);
306  }
307 
308  // The ECAL Properties
310  showerparam(myCalorimeter_->ecalProperties(onEcal),
311  myCalorimeter_->hcalProperties(onHcal),
312  myCalorimeter_->layer1Properties(onLayer1),
313  myCalorimeter_->layer2Properties(onLayer2),
316  RCFactor_,
317  RTFactor_);
318 
319  // Photons : create an e+e- pair
320  if ( myTrack.type() == 22 ) {
321 
322  // Depth for the first e+e- pair creation (in X0)
323  X0depth = -log(random->flatShoot()) * (9./7.);
324 
325  // Initialization
326  double eMass = 0.000510998902;
327  double xe=0;
328  double xm=eMass/myPart.e();
329  double weight = 0.;
330 
331  // Generate electron energy between emass and eGamma-emass
332  do {
333  xe = random->flatShoot()*(1.-2.*xm) + xm;
334  weight = 1. - 4./3.*xe*(1.-xe);
335  } while ( weight < random->flatShoot() );
336 
337  // Protection agains infinite loop in Famos Shower
338  if ( myPart.e()*xe < 0.055 || myPart.e()*(1.-xe) < 0.055 ) {
339 
340  if ( myPart.e() > 0.055 ) thePart.push_back(&myPart);
341 
342  } else {
343 
344  myElec = (myPart) * xe;
345  myPosi = (myPart) * (1.-xe);
348  thePart.push_back(&myElec);
349  thePart.push_back(&myPosi);
350  }
351  // Electrons
352  } else {
353 
354  X0depth = 0.;
355  if ( myPart.e() > 0.055 ) thePart.push_back(&myPart);
356 
357  }
358 
359  // After the different protections, this shouldn't happen.
360  if(thePart.empty())
361  {
362  if(myPreshower==nullptr) return;
363  delete myPreshower;
364  return;
365  }
366 
367  // find the most energetic particle
368  double maxEnergy=-1.;
369  for(unsigned ip=0;ip < thePart.size();++ip)
370  if(thePart[ip]->e() > maxEnergy) maxEnergy = thePart[ip]->e();
371 
372  // Initialize the Grid in ECAL
373  int size = gridSize_;
374  if(maxEnergy>100) size=11;
375  // if ( maxEnergy < threshold5x5 ) size = 5;
376  // if ( maxEnergy < threshold3x3 ) size = 3;
377 
378 
379  EMShower theShower(random,aGammaGenerator,&showerparam,&thePart,nullptr,nullptr,bFixedLength_);
380 
381 
382  double maxShower = theShower.getMaximumOfShower();
383  if (maxShower > 20.) maxShower = 2.; // simple pivot-searching protection
384 
385  double depth((X0depth + maxShower) *
387  XYZPoint meanShower = ecalentrance + myPart.Vect().Unit()*depth;
388 
389  // if(onEcal!=1) return ;
390 
391  // The closest crystal
392  DetId pivot(myCalorimeter_->getClosestCell(meanShower, true, onEcal==1));
393 
394  if(pivot.subdetId() == 0) { // further protection against avbsence of pivot
395  edm::LogWarning("CalorimetryManager") << "Pivot for egamma e = " << myTrack.hcalEntrance().e() << " is not found at depth " << depth << " and meanShower coordinates = " << meanShower << std::endl;
396  if(myPreshower) delete myPreshower;
397  return;
398  }
399 
400  EcalHitMaker myGrid(myCalorimeter_,ecalentrance,pivot,onEcal,size,0,random);
401  // ^^^^
402  // for EM showers
405 
406  //maximumdepth dependence of the radiusfactorbehindpreshower
407  //First tuning: Shilpi Jain (Mar-Apr 2010); changed after tuning - Feb-July - Shilpi Jain
408  /* **************
409  myGrid.setRadiusFactor(radiusFactor_);
410  if(onLayer1 || onLayer2)
411  {
412  float b = radiusPreshowerCorrections_[0];
413  float a = radiusFactor_*( 1.+radiusPreshowerCorrections_[1]*radiusPreshowerCorrections_[0] );
414  float maxdepth = X0depth+theShower.getMaximumOfShower();
415  float newRadiusFactor = radiusFactor_;
416  if(myPart.e()<=250.)
417  {
418  newRadiusFactor = a/(1.+b*maxdepth);
419  }
420  myGrid.setRadiusFactor(newRadiusFactor);
421  }
422  else // otherwise use the normal radius factor
423  {
424  myGrid.setRadiusFactor(radiusFactor_);
425  }
426  ************** */
427  if(myTrack.onEcal() == 2) // if on EE
428  {
429  if( (onLayer1 || onLayer2) && myPart.e()<=250.)
430  {
431  double maxdepth = X0depth+theShower.getMaximumOfShower();
432  double newRadiusFactor = radiusFactorEE_ * aTerm/(1.+bTerm*maxdepth);
433  myGrid.setRadiusFactor(newRadiusFactor);
434  }
435  else // otherwise use the normal radius factor
436  {
438  }
439  }//if(myTrack.onEcal() == 2)
440  else // else if on EB
441  {
443  }
444  //(end of) changed after tuning - Feb-July - Shilpi Jain
445 
447 
448  // The shower simulation
449  myGrid.setTrackParameters(myPart.Vect().Unit(),X0depth,myTrack);
450 
451  // std::cout << " PS ECAL GAP HCAL X0 " << myGrid.ps1TotalX0()+myGrid.ps2TotalX0() << " " << myGrid.ecalTotalX0();
452  // std::cout << " " << myGrid.ecalHcalGapTotalX0() << " " << myGrid.hcalTotalX0() << std::endl;
453  // std::cout << " PS ECAL GAP HCAL L0 " << myGrid.ps1TotalL0()+myGrid.ps2TotalL0() << " " << myGrid.ecalTotalL0();
454  // std::cout << " " << myGrid.ecalHcalGapTotalL0() << " " << myGrid.hcalTotalL0() << std::endl;
455  // std::cout << "ECAL-HCAL " << myTrack.momentum().eta() << " " << myGrid.ecalHcalGapTotalL0() << std::endl;
456  //
457  // std::cout << " Grid created " << std::endl;
458  if(myPreshower) theShower.setPreshower(myPreshower);
459 
460  HcalHitMaker myHcalHitMaker(myGrid,(unsigned)0);
461 
462  theShower.setGrid(&myGrid);
463  theShower.setHcal(&myHcalHitMaker);
464  theShower.compute();
465  //myHistos->fill("h502", myPart->eta(),myGrid.totalX0());
466 
467  // calculate the total simulated energy for this particle
468  float simE = 0;
469  for( const auto& mapIterator : myGrid.getHits() ) {
470  simE += mapIterator.second;
471  }
472 
473  auto scale = ecalCorrection ? ecalCorrection->getScale( myTrack.ecalEntrance().e(),
474  std::abs( myTrack.ecalEntrance().eta() ), simE ) : 1.;
475 
476  // Save the hits !
477  updateECAL( myGrid.getHits(), onEcal,myTrack.id(), scale );
478 
479  // Now fill the HCAL hits
480  updateHCAL(myHcalHitMaker.getHits(),myTrack.id());
481 
482  // delete the preshower
483  if(myPreshower!=nullptr) {
484  updatePreshower(myPreshower->getHits(),myTrack.id());
485  delete myPreshower;
486  // std::cout << " Deleting myPreshower " << std::endl;
487  }
488 
489 }
490 
493 {
494  int hit;
495  int pid = abs(myTrack.type());
496  if (debug_) {
497  LogInfo("FastCalorimetry") << " reconstructHCAL " << myTrack << std::endl;
498  }
499 
500  // FSimTrack myTrack = mySimEvent.track(fsimi);
501 
502  // int pid=abs(myTrack.type());
503  // std::cout << "reconstructHCAL " << std::endl;
504 
505  XYZTLorentzVector trackPosition;
506  if (myTrack.onHcal()) {
507  trackPosition=myTrack.hcalEntrance().vertex();
508  hit = myTrack.onHcal()-1;
509  } else {
510  trackPosition=myTrack.vfcalEntrance().vertex();
511  hit = 2;
512  }
513 
514  double pathEta = trackPosition.eta();
515  double pathPhi = trackPosition.phi();
516  // double pathTheta = trackPosition.theta();
517 
518  double EGen = myTrack.hcalEntrance().e();
519  double emeas = 0.;
520  //double emeas = -0.0001;
521 
522  if(pid == 13) {
523  // std::cout << " We should not be here " << std::endl;
524  emeas = myHDResponse_->responseHCAL(0, EGen, pathEta, 2, random); // 2=muon
525  if(debug_)
526  LogInfo("FastCalorimetry") << "CalorimetryManager::reconstructHCAL - MUON !!!" << std::endl;
527  }
528  else if( pid == 22 || pid == 11) {
529  emeas = myHDResponse_->responseHCAL(0, EGen, pathEta, 0, random); // last par. = 0 = e/gamma
530  // cout << "CalorimetryManager::reconstructHCAL - e/gamma !!!" << std::endl;
531  if(debug_)
532  LogInfo("FastCalorimetry") << "CalorimetryManager::reconstructHCAL - e/gamma !!!" << std::endl;
533  }
534  else {
535  emeas = myHDResponse_->getHCALEnergyResponse(EGen, hit, random);
536  }
537 
538  if(debug_)
539  LogInfo("FastCalorimetry") << "CalorimetryManager::reconstructHCAL - on-calo "
540  << " eta = " << pathEta
541  << " phi = " << pathPhi
542  << " Egen = " << EGen
543  << " Emeas = " << emeas << std::endl;
544 
545  if(emeas > 0.) {
546  DetId cell = myCalorimeter_->getClosestCell(trackPosition.Vect(),false,false);
547  double tof = (((HcalGeometry*)(myCalorimeter_->getHcalGeometry()))->getPosition(cell).mag())/29.98;//speed of light
548  CaloHitID current_id(cell.rawId(),tof,myTrack.id());
549  std::map<CaloHitID,float> hitMap;
550  hitMap[current_id] = emeas;
551  updateHCAL(hitMap,myTrack.id());
552  }
553 }
554 
556  // const edm::ParameterSet& fastCalo){
557 
558  // TimeMe t(" FASTEnergyReconstructor::HDShower");
559  const XYZTLorentzVector& moment = myTrack.momentum();
560 
561  if(debug_)
562  LogInfo("FastCalorimetry")
563  << "CalorimetryManager::HDShowerSimulation - track param."
564  << std::endl
565  << " eta = " << moment.eta() << std::endl
566  << " phi = " << moment.phi() << std::endl
567  << " et = " << moment.Et() << std::endl
568  << " e = " << myTrack.hcalEntrance().e() << std::endl;
569 
570  if (debug_) {
571  LogInfo("FastCalorimetry") << " HDShowerSimulation " << myTrack << std::endl;
572  }
573 
574 
575  int hit;
576  // int pid = abs(myTrack.type());
577 
578  XYZTLorentzVector trackPosition;
579  if ( myTrack.onEcal() ) {
580  trackPosition=myTrack.ecalEntrance().vertex();
581  hit = myTrack.onEcal()-1; //
582  myPart = myTrack.ecalEntrance();
583  } else if ( myTrack.onVFcal()) {
584  trackPosition=myTrack.vfcalEntrance().vertex();
585  hit = 2;
586  myPart = myTrack.vfcalEntrance();
587  }
588  else
589  {
590  LogInfo("FastCalorimetry") << " The particle is not in the acceptance " << std::endl;
591  return;
592  }
593 
594  // int onHCAL = hit + 1; - specially for myCalorimeter->hcalProperties(onHCAL)
595  // (below) to get VFcal properties ...
596  int onHCAL = hit + 1;
597  int onECAL = myTrack.onEcal();
598 
599  double pathEta = trackPosition.eta();
600  double pathPhi = trackPosition.phi();
601  // double pathTheta = trackPosition.theta();
602 
603  double eint = moment.e();
604  double eGen = myTrack.hcalEntrance().e();
605 
606  double emeas = 0.;
607  double pmip= myHDResponse_->getMIPfraction(eGen, pathEta);
608  // std::cout << " CalorimetryManager onHcal " << pmip << std::endl;
609 
610 
611  //double emeas = -0.000001;
612 
613  //===========================================================================
614  if(eGen > 0.) {
615 
616  // ECAL and HCAL properties to get
618  theHDShowerparam(myCalorimeter_->ecalProperties(onECAL),
621 
622  //Making ECAL Grid (and segments calculation)
623  XYZPoint caloentrance;
624  XYZVector direction;
625  if(myTrack.onEcal())
626  {
627  caloentrance = myTrack.ecalEntrance().vertex().Vect();
628  direction = myTrack.ecalEntrance().Vect().Unit();
629  }
630  else if(myTrack.onHcal())
631  {
632  caloentrance = myTrack.hcalEntrance().vertex().Vect();
633  direction = myTrack.hcalEntrance().Vect().Unit();
634  }
635  else
636  {
637  caloentrance = myTrack.vfcalEntrance().vertex().Vect();
638  direction = myTrack.vfcalEntrance().Vect().Unit();
639  }
640 
641  if(debug_)
642  LogInfo("FastCalorimetry")
643  << "CalorimetryManager::HDShowerSimulation - on-calo 1 "
644  << std::endl
645  << " onEcal = " << myTrack.onEcal() << std::endl
646  << " onHcal = " << myTrack.onHcal() << std::endl
647  << " onVFcal = " << myTrack.onVFcal() << std::endl
648  << " position = " << caloentrance << std::endl;
649 
650 
651  DetId pivot;
652  if(myTrack.onEcal())
653  {
654  pivot=myCalorimeter_->getClosestCell(caloentrance,
655  true, myTrack.onEcal()==1);
656  }
657  else if(myTrack.onHcal())
658  {
659  // std::cout << " CalorimetryManager onHcal " << myTrack.onHcal() << " caloentrance" << caloentrance << std::endl;
660  pivot=myCalorimeter_->getClosestCell(caloentrance,
661  false, false);
662  }
663 
664  EcalHitMaker myGrid(myCalorimeter_,caloentrance,pivot,
665  pivot.null()? 0 : myTrack.onEcal(),hdGridSize_,1,
666  random);
667  // 1=HAD shower
668 
669  myGrid.setTrackParameters(direction,0,myTrack);
670  // Build the FAMOS HCAL
671  HcalHitMaker myHcalHitMaker(myGrid,(unsigned)1);
672 
673  // Shower simulation
674  bool status = false;
675  int mip = 2;
676  // Use HFShower for HF
677  if ( !myTrack.onEcal() && !myTrack.onHcal() ) {
678  // std::cout << "CalorimetryManager::HDShowerSimulation(): track entrance = "
679  // << myTrack.vfcalEntrance().vertex().X() << " "
680  // << myTrack.vfcalEntrance().vertex().Y() << " "
681  // << myTrack.vfcalEntrance().vertex().Z() << " "
682  // << " , Energy (Gen/Scale) = " << eGen << " " << e << std::endl;
683 
684  // Warning : We give here the particle energy with the response
685  // but without the resolution/gaussian smearing
686  // For HF, the resolution is due to the PE statistic
687 
688  if(useShowerLibrary) {
690  status = true;
691  } else {
692  HFShower theShower(random,
693  &theHDShowerparam,
694  &myGrid,
695  &myHcalHitMaker,
696  onECAL,
697  eGen);
698  // eGen);
699  // e); // PV Warning : temporarly set the energy to the generated E
700 
701  status = theShower.compute();
702  }
703  } else {
704  if(hdSimMethod_ == 0) {
705  HDShower theShower(random,
706  &theHDShowerparam,
707  &myGrid,
708  &myHcalHitMaker,
709  onECAL,
710  eGen,
711  pmip);
712  status = theShower.compute();
713  mip = theShower.getmip();
714  }
715  else if (hdSimMethod_ == 1) {
716  HDRShower theShower(random,
717  &theHDShowerparam,
718  &myGrid,
719  &myHcalHitMaker,
720  onECAL,
721  eGen);
722  status = theShower.computeShower();
723  mip = 2;
724  }
725  else if (hdSimMethod_ == 2 ) {
726  // std::cout << "Using GflashHadronShowerProfile hdSimMethod_ == 2" << std::endl;
727 
728  //dynamically loading a corresponding profile by the particle type
729  int particleType = myTrack.type();
731  if(particleType == -2212) theProfile = theAntiProtonProfile;
732  else if(particleType == 2212) theProfile = theProtonProfile;
733 
734  //input variables for GflashHadronShowerProfile
735  int showerType = 99 + myTrack.onEcal();
736  double globalTime = 150.0; // a temporary reference hit time in nanosecond
737  float charge = (float)(myTrack.charge());
738  Gflash3Vector gfpos(trackPosition.X(),trackPosition.Y(),trackPosition.Z());
739  Gflash3Vector gfmom(moment.X(),moment.Y(),moment.Z());
740 
741  theProfile->initialize(showerType,eGen,globalTime,charge,gfpos,gfmom);
744 
745  //make hits
746  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
747  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
748  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
749 
751 
752  for( ; spotIter != spotIterEnd; spotIter++){
753 
754  double pathLength = theProfile->getGflashShowino()->getPathLengthAtShower()
755  + (30*100/eGen)*(spotIter->getTime() - globalTime);
756 
757  double currentDepth = std::max(0.0,pathLength - theProfile->getGflashShowino()->getPathLengthOnEcal());
758 
759  //find the the showino position at the currentDepth
760  GflashTrajectoryPoint trajectoryPoint;
761  theProfile->getGflashShowino()->getHelix()->getGflashTrajectoryPoint(trajectoryPoint,pathLength);
762  Gflash3Vector positionAtCurrentDepth = trajectoryPoint.getPosition();
763  //find radial distrance
764  Gflash3Vector lateralDisplacement = positionAtCurrentDepth - spotIter->getPosition()/CLHEP::cm;
765  double rShower = lateralDisplacement.r();
766  double azimuthalAngle = lateralDisplacement.phi();
767 
768  whichCalor = Gflash::getCalorimeterNumber(positionAtCurrentDepth);
769 
770  if(whichCalor==Gflash::kESPM || whichCalor==Gflash::kENCA) {
771  bool statusPad = myGrid.getPads(currentDepth,true);
772  if(!statusPad) continue;
773  myGrid.setSpotEnergy(1.2*spotIter->getEnergy()/CLHEP::GeV);
774  myGrid.addHit(rShower/Gflash::intLength[Gflash::kESPM],azimuthalAngle,0);
775  }
776  else if(whichCalor==Gflash::kHB || whichCalor==Gflash::kHE) {
777  bool setHDdepth = myHcalHitMaker.setDepth(currentDepth,true);
778  if(!setHDdepth) continue;
779  myHcalHitMaker.setSpotEnergy(1.4*spotIter->getEnergy()/CLHEP::GeV);
780  myHcalHitMaker.addHit(rShower/Gflash::intLength[Gflash::kHB],azimuthalAngle,0);
781  }
782  }
783  status = true;
784  }
785  else {
786  edm::LogInfo("FastSimulationCalorimetry") << " SimMethod " << hdSimMethod_ <<" is NOT available ";
787  }
788  }
789 
790 
791  if(status) {
792 
793  // Here to switch between simple formulae and parameterized response
794  if(optionHDSim_ == 1) {
795  emeas = myHDResponse_->getHCALEnergyResponse(eGen, hit, random);
796  }
797  else { // optionHDsim == 2
798  emeas = myHDResponse_->responseHCAL(mip, eGen, pathEta, 1, random); // 1=hadron
799  }
800 
801  double correction = emeas / eGen;
802 
803  // RespCorrP factors (ECAL and HCAL separately) calculation
804  respCorr(eint);
805 
806  if(debug_)
807  LogInfo("FastCalorimetry")
808  << "CalorimetryManager::HDShowerSimulation - on-calo 2" << std::endl
809  << " eta = " << pathEta << std::endl
810  << " phi = " << pathPhi << std::endl
811  << " Egen = " << eGen << std::endl
812  << " Emeas = " << emeas << std::endl
813  << " corr = " << correction << std::endl
814  << " mip = " << mip << std::endl;
815 
816  if(myTrack.onEcal() > 0) {
817  // Save ECAL hits
818  updateECAL(myGrid.getHits(),onECAL,myTrack.id(),correction*ecorr);
819  }
820 
821  // Save HCAL hits
822  if(myTrack.onVFcal() && useShowerLibrary) {
823  myHDResponse_->correctHF(eGen,abs(myTrack.type()));
825  }
826  else
827  updateHCAL(myHcalHitMaker.getHits(),myTrack.id(),correction*hcorr);
828 
829  }
830  else { // shower simulation failed
831  // std::cout << " Shower simulation failed " << trackPosition.Vect() << std::endl;
832  // std::cout << " The FSimTrack " << myTrack << std::endl;
833  // std::cout << " HF entrance on VFcal" << myTrack.onVFcal() << std::endl;
834  // std::cout << " trackPosition.eta() " << trackPosition.eta() << std::endl;
835  if(myTrack.onHcal() || myTrack.onVFcal())
836  {
837  DetId cell = myCalorimeter_->getClosestCell(trackPosition.Vect(),false,false);
838  double tof = (((HcalGeometry*)(myCalorimeter_->getHcalGeometry()))->getPosition(cell).mag())/29.98;//speed of light
839  CaloHitID current_id(cell.rawId(),tof,myTrack.id());
840  std::map<CaloHitID,float> hitMap;
841  hitMap[current_id] = emeas;
842  updateHCAL(hitMap,myTrack.id());
843  if(debug_)
844  LogInfo("FastCalorimetry") << " HCAL simple cell "
845  << cell.rawId() << " added E = "
846  << emeas << std::endl;
847  }
848  }
849 
850  } // e > 0. ...
851 
852  if(debug_)
853  LogInfo("FastCalorimetry") << std::endl << " FASTEnergyReconstructor::HDShowerSimulation finished "
854  << std::endl;
855 }
856 
857 
859 {
860  // TimeMe t(" FASTEnergyReconstructor::HDShower");
861  XYZTLorentzVector moment = myTrack.momentum();
862 
863  // Backward compatibility behaviour
864  if(!theMuonHcalEffects)
865  {
866  if(myTrack.onHcal() || myTrack.onVFcal() )
867  reconstructHCAL(myTrack, random);
868 
869  return;
870  }
871 
872  if(debug_)
873  LogInfo("FastCalorimetry") << "CalorimetryManager::MuonMipSimulation - track param."
874  << std::endl
875  << " eta = " << moment.eta() << std::endl
876  << " phi = " << moment.phi() << std::endl
877  << " et = " << moment.Et() << std::endl;
878 
879  // int hit;
880  // int pid = abs(myTrack.type());
881 
882  XYZTLorentzVector trackPosition;
883  if ( myTrack.onEcal() ) {
884  trackPosition=myTrack.ecalEntrance().vertex();
885  // hit = myTrack.onEcal()-1; //
886  myPart = myTrack.ecalEntrance();
887  } else if ( myTrack.onVFcal()) {
888  trackPosition=myTrack.vfcalEntrance().vertex();
889  // hit = 2;
890  myPart = myTrack.vfcalEntrance();
891  }
892  else
893  {
894  LogInfo("FastCalorimetry") << " The particle is not in the acceptance " << std::endl;
895  return;
896  }
897 
898  // int onHCAL = hit + 1; - specially for myCalorimeter->hcalProperties(onHCAL)
899  // (below) to get VFcal properties ...
900  // not needed ?
901  // int onHCAL = hit + 1;
902  int onECAL = myTrack.onEcal();
903 
904  // double pathEta = trackPosition.eta();
905  // double pathPhi = trackPosition.phi();
906  // double pathTheta = trackPosition.theta();
907 
908  //===========================================================================
909 
910  // ECAL and HCAL properties to get
911 
912  //Making ECAL Grid (and segments calculation)
913  XYZPoint caloentrance;
914  XYZVector direction;
915  if(myTrack.onEcal())
916  {
917  caloentrance = myTrack.ecalEntrance().vertex().Vect();
918  direction = myTrack.ecalEntrance().Vect().Unit();
919  }
920  else if(myTrack.onHcal())
921  {
922  caloentrance = myTrack.hcalEntrance().vertex().Vect();
923  direction = myTrack.hcalEntrance().Vect().Unit();
924  }
925  else
926  {
927  caloentrance = myTrack.vfcalEntrance().vertex().Vect();
928  direction = myTrack.vfcalEntrance().Vect().Unit();
929  }
930 
931  DetId pivot;
932  if(myTrack.onEcal())
933  {
934  pivot=myCalorimeter_->getClosestCell(caloentrance,
935  true, myTrack.onEcal()==1);
936  }
937  else if(myTrack.onHcal())
938  {
939  // std::cout << " CalorimetryManager onHcal " << myTrack.onHcal() << " caloentrance" << caloentrance << std::endl;
940  pivot=myCalorimeter_->getClosestCell(caloentrance,
941  false, false);
942  }
943 
944  EcalHitMaker myGrid(myCalorimeter_,caloentrance,pivot,
945  pivot.null()? 0 : myTrack.onEcal(),hdGridSize_,0,
946  random);
947  // 0 =EM shower -> Unit = X0
948 
949  myGrid.setTrackParameters(direction,0,myTrack);
950 
951  // Now get the path in the Preshower, ECAL and HCAL along a straight line extrapolation
952  // but only those in the ECAL are used
953 
954  const std::vector<CaloSegment>& segments=myGrid.getSegments();
955  unsigned nsegments=segments.size();
956 
957  int ifirstHcal=-1;
958  int ilastEcal=-1;
959 
960  EnergyLossSimulator* energyLossECAL = (theMuonEcalEffects) ?
962  // // Muon brem in ECAL
963  // MuonBremsstrahlungSimulator* muonBremECAL = 0;
964  // if (theMuonEcalEffects) muonBremECAL = theMuonEcalEffects->muonBremsstrahlungSimulator();
965 
966  for(unsigned iseg=0;iseg<nsegments&&ifirstHcal<0;++iseg)
967  {
968 
969  // in the ECAL, there are two types of segments: PbWO4 and GAP
970  float segmentSizeinX0=segments[iseg].X0length();
971 
972  // Martijn - insert your computations here
973  float energy=0.0;
974  if (segmentSizeinX0>0.001 && segments[iseg].material()==CaloSegment::PbWO4 ) {
975  // The energy loss simulator
976  float charge = (float)(myTrack.charge());
977  ParticlePropagator theMuon(moment,trackPosition,charge,nullptr);
978  theMuon.setID(-(int)charge*13);
979  if ( energyLossECAL ) {
980  energyLossECAL->updateState(theMuon, segmentSizeinX0, random);
981  energy = energyLossECAL->deltaMom().E();
982  moment -= energyLossECAL->deltaMom();
983  }
984  }
985  // that's all for ECAL, Florian
986  // Save the hit only if it is a crystal
987  if(segments[iseg].material()==CaloSegment::PbWO4)
988  {
989  myGrid.getPads(segments[iseg].sX0Entrance()+segmentSizeinX0*0.5);
990  myGrid.setSpotEnergy(energy);
991  myGrid.addHit(0.,0.);
992  ilastEcal=iseg;
993  }
994  // Check for end of loop:
995  if(segments[iseg].material()==CaloSegment::HCAL)
996  {
997  ifirstHcal=iseg;
998  }
999  }
1000 
1001 
1002  // Build the FAMOS HCAL
1003  HcalHitMaker myHcalHitMaker(myGrid,(unsigned)2);
1004  // float mipenergy=0.1;
1005  // Create the helix with the stepping helix propagator
1006  // to add a hit, just do
1007  // myHcalHitMaker.setSpotEnergy(mipenergy);
1008  // math::XYZVector hcalEntrance;
1009  // if(ifirstHcal>=0) hcalEntrance=segments[ifirstHcal].entrance();
1010  // myHcalHitMaker.addHit(hcalEntrance);
1014  int ilastHcal=-1;
1015  float mipenergy=0.0;
1016 
1017  EnergyLossSimulator* energyLossHCAL = (theMuonHcalEffects) ?
1019  // // Muon Brem effect
1020  // MuonBremsstrahlungSimulator* muonBremHCAL = 0;
1021  // if (theMuonHcalEffects) muonBremHCAL = theMuonHcalEffects->muonBremsstrahlungSimulator();
1022 
1023  if(ifirstHcal>0 && energyLossHCAL){
1024  for(unsigned iseg=ifirstHcal;iseg<nsegments;++iseg)
1025  {
1026  float segmentSizeinX0=segments[iseg].X0length();
1027  if(segments[iseg].material()==CaloSegment::HCAL) {
1028  ilastHcal=iseg;
1029  if (segmentSizeinX0>0.001) {
1030  // The energy loss simulator
1031  float charge = (float)(myTrack.charge());
1032  ParticlePropagator theMuon(moment,trackPosition,charge,nullptr);
1033  theMuon.setID(-(int)charge*13);
1034  energyLossHCAL->updateState(theMuon, segmentSizeinX0, random);
1035  mipenergy = energyLossHCAL->deltaMom().E();
1036  moment -= energyLossHCAL->deltaMom();
1037  myHcalHitMaker.setSpotEnergy(mipenergy);
1038  myHcalHitMaker.addHit(segments[iseg].entrance());
1039  }
1040  }
1041  }
1042  }
1043 
1048  //
1049 
1050 
1051 
1052  // Copy the muon SimTrack (Only for Energy loss)
1053  FSimTrack muonTrack(myTrack);
1054  if(energyLossHCAL && ilastHcal>=0) {
1055  math::XYZVector hcalExit=segments[ilastHcal].exit();
1056  muonTrack.setTkPosition(hcalExit);
1057  muonTrack.setTkMomentum(moment);
1058  } else if(energyLossECAL && ilastEcal>=0) {
1059  math::XYZVector ecalExit=segments[ilastEcal].exit();
1060  muonTrack.setTkPosition(ecalExit);
1061  muonTrack.setTkMomentum(moment);
1062  } // else just leave tracker surface position and momentum...
1063 
1064  muonSimTracks.push_back(muonTrack);
1065 
1066 
1067  // no need to change below this line
1068  std::map<CaloHitID,float>::const_iterator mapitr;
1069  std::map<CaloHitID,float>::const_iterator endmapitr;
1070  if(myTrack.onEcal() > 0) {
1071  // Save ECAL hits
1072  updateECAL(myGrid.getHits(),onECAL,myTrack.id());
1073  }
1074 
1075  // Save HCAL hits
1076  updateHCAL(myHcalHitMaker.getHits(),myTrack.id());
1077 
1078  if(debug_)
1079  LogInfo("FastCalorimetry") << std::endl << " FASTEnergyReconstructor::MipShowerSimulation finished "
1080  << std::endl;
1081 }
1082 
1083 
1085 
1086  edm::ParameterSet ECALparameters = fastCalo.getParameter<edm::ParameterSet>("ECAL");
1087 
1088  evtsToDebug_ = fastCalo.getUntrackedParameter<std::vector<unsigned int> >("EvtsToDebug",std::vector<unsigned>());
1089  debug_ = fastCalo.getUntrackedParameter<bool>("Debug");
1090 
1091  bFixedLength_ = ECALparameters.getParameter<bool>("bFixedLength");
1092  // std::cout << "bFixedLength_ = " << bFixedLength_ << std::endl;
1093 
1094  gridSize_ = ECALparameters.getParameter<int>("GridSize");
1095  spotFraction_ = ECALparameters.getParameter<double>("SpotFraction");
1096  pulledPadSurvivalProbability_ = ECALparameters.getParameter<double>("FrontLeakageProbability");
1097  crackPadSurvivalProbability_ = ECALparameters.getParameter<double>("GapLossProbability");
1098  theCoreIntervals_ = ECALparameters.getParameter<std::vector<double> >("CoreIntervals");
1099  theTailIntervals_ = ECALparameters.getParameter<std::vector<double> >("TailIntervals");
1100 
1101  RCFactor_ = ECALparameters.getParameter<double>("RCFactor");
1102  RTFactor_ = ECALparameters.getParameter<double>("RTFactor");
1103  //changed after tuning - Feb-July - Shilpi Jain
1104  // radiusFactor_ = ECALparameters.getParameter<double>("RadiusFactor");
1105  radiusFactorEE_ = ECALparameters.getParameter<double>("RadiusFactorEE");
1106  radiusFactorEB_ = ECALparameters.getParameter<double>("RadiusFactorEB");
1107  //(end of) changed after tuning - Feb-July - Shilpi Jain
1108  radiusPreshowerCorrections_ = ECALparameters.getParameter<std::vector<double> >("RadiusPreshowerCorrections");
1109  aTerm = 1.+radiusPreshowerCorrections_[1]*radiusPreshowerCorrections_[0];
1110  bTerm = radiusPreshowerCorrections_[0];
1111  mipValues_ = ECALparameters.getParameter<std::vector<double> >("MipsinGeV");
1112  simulatePreshower_ = ECALparameters.getParameter<bool>("SimulatePreshower");
1113 
1114  if(gridSize_ <1) gridSize_= 7;
1115  if(pulledPadSurvivalProbability_ <0. || pulledPadSurvivalProbability_>1 ) pulledPadSurvivalProbability_= 1.;
1116  if(crackPadSurvivalProbability_ <0. || crackPadSurvivalProbability_>1 ) crackPadSurvivalProbability_= 0.9;
1117 
1118  LogInfo("FastCalorimetry") << " Fast ECAL simulation parameters " << std::endl;
1119  LogInfo("FastCalorimetry") << " =============================== " << std::endl;
1120  if(simulatePreshower_)
1121  LogInfo("FastCalorimetry") << " The preshower is present " << std::endl;
1122  else
1123  LogInfo("FastCalorimetry") << " The preshower is NOT present " << std::endl;
1124  LogInfo("FastCalorimetry") << " Grid Size : " << gridSize_ << std::endl;
1125  if(spotFraction_>0.)
1126  LogInfo("FastCalorimetry") << " Spot Fraction : " << spotFraction_ << std::endl;
1127  else
1128  {
1129  LogInfo("FastCalorimetry") << " Core of the shower " << std::endl;
1130  for(unsigned ir=0; ir < theCoreIntervals_.size()/2;++ir)
1131  {
1132  LogInfo("FastCalorimetry") << " r < " << theCoreIntervals_[ir*2] << " R_M : " << theCoreIntervals_[ir*2+1] << " ";
1133  }
1134  LogInfo("FastCalorimetry") << std::endl;
1135 
1136  LogInfo("FastCalorimetry") << " Tail of the shower " << std::endl;
1137  for(unsigned ir=0; ir < theTailIntervals_.size()/2;++ir)
1138  {
1139  LogInfo("FastCalorimetry") << " r < " << theTailIntervals_[ir*2] << " R_M : " << theTailIntervals_[ir*2+1] << " ";
1140  }
1141  //changed after tuning - Feb-July - Shilpi Jain
1142  // LogInfo("FastCalorimetry") << "Radius correction factor " << radiusFactor_ << std::endl;
1143  LogInfo("FastCalorimetry") << "Radius correction factors: EB & EE " << radiusFactorEB_ << " : "<< radiusFactorEE_ << std::endl;
1144  //(end of) changed after tuning - Feb-July - Shilpi Jain
1145  LogInfo("FastCalorimetry") << std::endl;
1146  if(mipValues_.size()>2) {
1147  LogInfo("FastCalorimetry") << "Improper number of parameters for the preshower ; using 95keV" << std::endl;
1148  mipValues_.clear();
1149  mipValues_.resize(2,0.000095);
1150  }
1151  }
1152 
1153  LogInfo("FastCalorimetry") << " FrontLeakageProbability : " << pulledPadSurvivalProbability_ << std::endl;
1154  LogInfo("FastCalorimetry") << " GapLossProbability : " << crackPadSurvivalProbability_ << std::endl;
1155 
1156 
1157  // RespCorrP: p (momentum), ECAL and HCAL corrections = f(p)
1158  edm::ParameterSet CalorimeterParam = fastCalo.getParameter<edm::ParameterSet>("CalorimeterProperties");
1159 
1160  rsp = CalorimeterParam.getParameter<std::vector<double> >("RespCorrP");
1161  LogInfo("FastCalorimetry") << " RespCorrP (rsp) size " << rsp.size() << std::endl;
1162 
1163  if( rsp.size()%3 !=0 ) {
1164  LogInfo("FastCalorimetry")
1165  << " RespCorrP size is wrong -> no corrections applied !!!"
1166  << std::endl;
1167 
1168  p_knots.push_back(14000.);
1169  k_e.push_back (1.);
1170  k_h.push_back (1.);
1171  }
1172  else {
1173  for(unsigned i = 0; i < rsp.size(); i += 3) {
1174  LogInfo("FastCalorimetry") << "i = " << i/3 << " p = " << rsp [i]
1175  << " k_e(p) = " << rsp[i+1]
1176  << " k_e(p) = " << rsp[i+2] << std::endl;
1177 
1178  p_knots.push_back(rsp[i]);
1179  k_e.push_back (rsp[i+1]);
1180  k_h.push_back (rsp[i+2]);
1181  }
1182  }
1183 
1184 
1185  //FR
1186  edm::ParameterSet HCALparameters = fastCalo.getParameter<edm::ParameterSet>("HCAL");
1187  optionHDSim_ = HCALparameters.getParameter<int>("SimOption");
1188  hdGridSize_ = HCALparameters.getParameter<int>("GridSize");
1189  hdSimMethod_ = HCALparameters.getParameter<int>("SimMethod");
1190  //RF
1191 
1192  EcalDigitizer_ = ECALparameters.getUntrackedParameter<bool>("Digitizer",false);
1193  HcalDigitizer_ = HCALparameters.getUntrackedParameter<bool>("Digitizer",false);
1194  samplingHBHE_ = HCALparameters.getParameter< std::vector<double> >("samplingHBHE");
1195  samplingHF_ = HCALparameters.getParameter< std::vector<double> >("samplingHF");
1196  samplingHO_ = HCALparameters.getParameter< std::vector<double> >("samplingHO");
1197  ietaShiftHB_ = HCALparameters.getParameter< int >("ietaShiftHB");
1198  ietaShiftHE_ = HCALparameters.getParameter< int >("ietaShiftHE");
1199  ietaShiftHF_ = HCALparameters.getParameter< int >("ietaShiftHF");
1200  ietaShiftHO_ = HCALparameters.getParameter< int >("ietaShiftHO");
1201  timeShiftHB_ = HCALparameters.getParameter< std::vector<double> >("timeShiftHB");
1202  timeShiftHE_ = HCALparameters.getParameter< std::vector<double> >("timeShiftHE");
1203  timeShiftHF_ = HCALparameters.getParameter< std::vector<double> >("timeShiftHF");
1204  timeShiftHO_ = HCALparameters.getParameter< std::vector<double> >("timeShiftHO");
1205 
1206  // FastHFShowerLibrary
1207  edm::ParameterSet m_HS = fastCalo.getParameter<edm::ParameterSet>("HFShowerLibrary");
1208  useShowerLibrary = m_HS.getUntrackedParameter<bool>("useShowerLibrary",false);
1209  useCorrectionSL = m_HS.getUntrackedParameter<bool>("useCorrectionSL",false);
1210 }
1211 
1213 
1214  int sizeP = p_knots.size();
1215 
1216  if(sizeP <= 1) {
1217  ecorr = 1.;
1218  hcorr = 1.;
1219  }
1220  else {
1221  int ip = -1;
1222  for (int i = 0; i < sizeP; i++) {
1223  if (p < p_knots[i]) { ip = i; break;}
1224  }
1225  if (ip == 0) {
1226  ecorr = k_e[0];
1227  hcorr = k_h[0];
1228  }
1229  else {
1230  if(ip == -1) {
1231  ecorr = k_e[sizeP-1];
1232  hcorr = k_h[sizeP-1];
1233  }
1234  else {
1235  double x1 = p_knots[ip-1];
1236  double x2 = p_knots[ip];
1237  double y1 = k_e[ip-1];
1238  double y2 = k_e[ip];
1239 
1240  if(x1 == x2) {
1241  // std::cout << " equal p_knots values!!! " << std::endl;
1242  }
1243 
1244  ecorr = (y1 + (y2 - y1) * (p - x1)/(x2 - x1));
1245 
1246  y1 = k_h[ip-1];
1247  y2 = k_h[ip];
1248  hcorr = (y1 + (y2 - y1) * (p - x1)/(x2 - x1));
1249 
1250  }
1251  }
1252  }
1253 
1254  if(debug_)
1255  LogInfo("FastCalorimetry") << " p, ecorr, hcorr = " << p << " "
1256  << ecorr << " " << hcorr << std::endl;
1257 
1258 }
1259 
1260 void CalorimetryManager::updateECAL(const std::map<CaloHitID,float>& hitMap, int onEcal, int trackID, float corr)
1261 {
1262  std::map<CaloHitID,float>::const_iterator mapitr;
1263  std::map<CaloHitID,float>::const_iterator endmapitr=hitMap.end();
1264  if(onEcal==1) {
1265  EBMapping_.reserve(EBMapping_.size()+hitMap.size());
1266  endmapitr=hitMap.end();
1267  for(mapitr=hitMap.begin();mapitr!=endmapitr;++mapitr) {
1268  //correct energy
1269  float energy = mapitr->second;
1270  energy *= corr;
1271 
1272  //make finalized CaloHitID
1273  CaloHitID current_id(mapitr->first.unitID(),mapitr->first.timeSlice(),trackID);
1274 
1275  EBMapping_.push_back(std::pair<CaloHitID,float>(current_id,energy));
1276  }
1277  }
1278  else if(onEcal==2) {
1279  EEMapping_.reserve(EEMapping_.size()+hitMap.size());
1280  endmapitr=hitMap.end();
1281  for(mapitr=hitMap.begin();mapitr!=endmapitr;++mapitr) {
1282  //correct energy
1283  float energy = mapitr->second;
1284  energy *= corr;
1285 
1286  //make finalized CaloHitID
1287  CaloHitID current_id(mapitr->first.unitID(),mapitr->first.timeSlice(),trackID);
1288 
1289  EEMapping_.push_back(std::pair<CaloHitID,float>(current_id,energy));
1290  }
1291  }
1292 
1293 }
1294 
1295 void CalorimetryManager::updateHCAL(const std::map<CaloHitID,float>& hitMap, int trackID, float corr)
1296 {
1297  std::vector<double> hfcorrEm = myHDResponse_->getCorrHFem();
1298  std::vector<double> hfcorrHad = myHDResponse_->getCorrHFhad();
1299  std::map<CaloHitID,float>::const_iterator mapitr;
1300  std::map<CaloHitID,float>::const_iterator endmapitr=hitMap.end();
1301  HMapping_.reserve(HMapping_.size()+hitMap.size());
1302  for(mapitr=hitMap.begin(); mapitr!=endmapitr; ++mapitr) {
1303  //correct energy
1304  float energy = mapitr->second;
1305  energy *= corr;
1306 
1307  float time = mapitr->first.timeSlice();
1308  //put energy into uncalibrated state for digitizer && correct timing
1309  if(HcalDigitizer_){
1310  HcalDetId hdetid = HcalDetId(mapitr->first.unitID());
1311  if (hdetid.subdetId()== HcalBarrel){
1312  energy /= samplingHBHE_[hdetid.ietaAbs()-1]; //re-convert to GeV
1313  time = timeShiftHB_[hdetid.ietaAbs()-ietaShiftHB_];
1314  }
1315  else if (hdetid.subdetId()== HcalEndcap){
1316  energy /= samplingHBHE_[hdetid.ietaAbs()-1]; //re-convert to GeV
1317  time = timeShiftHE_[hdetid.ietaAbs()-ietaShiftHE_];
1318  }
1319  else if (hdetid.subdetId()== HcalForward){
1320  if(useShowerLibrary) {
1321  if(useCorrectionSL) {
1322  if(hdetid.depth()== 1) energy *= hfcorrEm[hdetid.ietaAbs()-ietaShiftHF_];
1323  if(hdetid.depth()== 2) energy *= hfcorrHad[hdetid.ietaAbs()-ietaShiftHF_];
1324  }
1325  } else {
1326  if(hdetid.depth()== 1) energy *= samplingHF_[0];
1327  if(hdetid.depth()== 2) energy *= samplingHF_[1];
1328  time = timeShiftHF_[hdetid.ietaAbs()-ietaShiftHF_];
1329  }
1330  }
1331  else if (hdetid.subdetId()== HcalOuter){
1332  energy /= samplingHO_[hdetid.ietaAbs()-1];
1333  time = timeShiftHO_[hdetid.ietaAbs()-ietaShiftHO_];
1334  }
1335  }
1336 
1337  //make finalized CaloHitID
1338  CaloHitID current_id(mapitr->first.unitID(),time,trackID);
1339  HMapping_.push_back(std::pair<CaloHitID,float>(current_id,energy));
1340  }
1341 }
1342 
1343 void CalorimetryManager::updatePreshower(const std::map<CaloHitID,float>& hitMap, int trackID, float corr)
1344 {
1345  std::map<CaloHitID,float>::const_iterator mapitr;
1346  std::map<CaloHitID,float>::const_iterator endmapitr=hitMap.end();
1347  ESMapping_.reserve(ESMapping_.size()+hitMap.size());
1348  for(mapitr=hitMap.begin();mapitr!=endmapitr;++mapitr) {
1349  //correct energy
1350  float energy = mapitr->second;
1351  energy *= corr;
1352 
1353  //make finalized CaloHitID
1354  CaloHitID current_id(mapitr->first.unitID(),mapitr->first.timeSlice(),trackID);
1355 
1356  ESMapping_.push_back(std::pair<CaloHitID,float>(current_id,energy));
1357  }
1358 }
1359 
1361 {
1362  c.reserve(c.size()+EBMapping_.size());
1363  for(unsigned i=0; i<EBMapping_.size(); i++) {
1364  c.push_back(PCaloHit(EBDetId::unhashIndex(EBMapping_[i].first.unitID()),EBMapping_[i].second,EBMapping_[i].first.timeSlice(),EBMapping_[i].first.trackID()));
1365  }
1366 }
1367 
1369 {
1370  c.reserve(c.size()+EEMapping_.size());
1371  for(unsigned i=0; i<EEMapping_.size(); i++) {
1372  c.push_back(PCaloHit(EEDetId::unhashIndex(EEMapping_[i].first.unitID()),EEMapping_[i].second,EEMapping_[i].first.timeSlice(),EEMapping_[i].first.trackID()));
1373  }
1374 }
1375 
1377 {
1378  c.reserve(c.size()+HMapping_.size());
1379  for(unsigned i=0; i<HMapping_.size(); i++) {
1380  c.push_back(PCaloHit(DetId(HMapping_[i].first.unitID()),HMapping_[i].second,HMapping_[i].first.timeSlice(),HMapping_[i].first.trackID()));
1381  }
1382 }
1383 
1384 
1386 {
1387  c.reserve(c.size()+ESMapping_.size());
1388  for(unsigned i=0; i<ESMapping_.size(); i++) {
1389  c.push_back(PCaloHit(ESMapping_[i].first.unitID(),ESMapping_[i].second,ESMapping_[i].first.timeSlice(),ESMapping_[i].first.trackID()));
1390  }
1391 }
1392 
1393 // The main danger in this method is to screw up to relationships between particles
1394 // So, the muon FSimTracks created by FSimEvent.cc are simply to be updated
1396 {
1397  unsigned size=muons.size();
1398  for(unsigned i=0; i<size;++i)
1399  {
1400  int id=muons[i].trackId();
1401  if(abs(muons[i].type())!=13) continue;
1402  // identify the corresponding muon in the local collection
1403 
1404  std::vector<FSimTrack>::const_iterator itcheck=find_if(muonSimTracks.begin(),muonSimTracks.end(),FSimTrackEqual(id));
1405  if(itcheck!=muonSimTracks.end())
1406  {
1407  muons[i].setTkPosition(itcheck->trackerSurfacePosition());
1408  muons[i].setTkMomentum(itcheck->trackerSurfaceMomentum());
1409  // std::cout << " Found the SimTrack " << std::endl;
1410  // std::cout << *itcheck << std::endl;
1411  // std::cout << "SimTrack Id "<< id << " " << muons[i] << " " << std::endl;
1412  }
1413  // else
1414  // {
1415  // std::cout << " Calorimetery Manager : this should really not happen " << std::endl;
1416  // std::cout << " Was looking for " << id << " " << muons[i] << std::endl;
1417  // for(unsigned i=0;i<muonSimTracks.size();++i)
1418  // std::cout << muonSimTracks[i] << std::endl;
1419  // }
1420  }
1421 
1422 }
1423 
size
Write out results.
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:41
T getUntrackedParameter(std::string const &, T const &) const
std::vector< double > rsp
std::vector< double > k_h
bool noEndVertex() const
no end vertex
RawParticle myElec
A few pointers to save time.
std::vector< std::pair< CaloHitID, float > > ESMapping_
float charge() const
charge
Definition: FSimTrack.h:47
void reconstructHCAL(const FSimTrack &myTrack, RandomEngineAndDistribution const *)
std::vector< PCaloHit > PCaloHitContainer
const double GeV
Definition: MathUtil.h:16
const ECALProperties * ecalProperties(int onEcal) const
ECAL properties.
Definition: Calorimeter.cc:76
double responseHCAL(int _mip, double energy, double eta, int partype, RandomEngineAndDistribution const *)
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:139
void updateHCAL(const std::map< CaloHitID, float > &hitMap, int trackID=0, float corr=1.0)
int onLayer2() const
Definition: FSimTrack.h:96
GflashPiKShowerProfile * thePiKProfile
double flatShoot(double xmin=0.0, double xmax=1.0) const
std::vector< double > timeShiftHO_
void updatePreshower(const std::map< CaloHitID, float > &hitMap, int trackID=0, float corr=1.0)
GflashTrajectory * getHelix()
Definition: GflashShowino.h:30
MaterialEffects * theMuonEcalEffects
MaterialEffects * theMuonHcalEffects
void recoHFShowerLibrary(const FSimTrack &myTrack)
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:190
void setCrackPadSurvivalProbability(double val)
Definition: EcalHitMaker.h:134
double getHCALEnergyResponse(double e, int hit, RandomEngineAndDistribution const *)
GflashHadronShowerProfile * theProfile
const RawParticle & layer1Entrance() const
The particle at Preshower Layer 1.
Definition: FSimTrack.h:127
bool compute()
Compute the shower longitudinal and lateral development.
Definition: HDShower.cc:447
bool exists(std::string const &parameterName) const
checks if a parameter exists
void correctHF(double e, int type)
void setPreshower(PreshowerHitMaker *const myPresh)
set the preshower address
Definition: EMShower.cc:704
void loadFromPreshower(edm::PCaloHitContainer &c) const
Definition: weight.py:1
void updateState(ParticlePropagator &myTrack, double radlen, RandomEngineAndDistribution const *)
Compute the material effect (calls the sub class)
void updateECAL(const std::map< CaloHitID, float > &hitMap, int onEcal, int trackID=0, float corr=1.0)
const std::map< CaloHitID, float > & getHits() override
const CaloSubdetectorGeometry * getHcalGeometry() const
Definition: Calorimeter.h:57
void MuonMipSimulation(const FSimTrack &myTrack, RandomEngineAndDistribution const *)
void readParameters(const edm::ParameterSet &fastCalo)
TRandom random
Definition: MVATrainer.cc:138
const PreshowerLayer1Properties * layer1Properties(int onLayer1) const
Preshower Layer1 properties.
Definition: Calorimeter.cc:103
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:99
#define nullptr
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
std::vector< double > timeShiftHF_
std::vector< double > p_knots
void compute()
Compute the shower longitudinal and lateral development.
Definition: EMShower.cc:273
std::vector< std::pair< CaloHitID, float > > EBMapping_
int onEcal() const
Definition: FSimTrack.h:101
std::vector< FSimTrack > muonSimTracks
double radLenIncm() const override
Radiation length in cm.
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
U second(std::pair< T, U > const &p)
void HDShowerSimulation(const FSimTrack &myTrack, RandomEngineAndDistribution const *)
Hadronic Shower Simulation.
const LandauFluctuationGenerator * aLandauGenerator
std::vector< double > samplingHF_
bool compute()
Compute the shower longitudinal and lateral development.
Definition: HFShower.cc:475
void setRadiusFactor(double r)
Definition: EcalHitMaker.h:130
bool addHit(double r, double phi, unsigned layer=0) override
add the hit in the HCAL in local coordinates
Definition: HcalHitMaker.cc:28
int depth() const
get the tower depth
Definition: HcalDetId.cc:108
void reconstruct(RandomEngineAndDistribution const *)
CalorimeterNumber getCalorimeterNumber(const Gflash3Vector &position)
GflashAntiProtonShowerProfile * theAntiProtonProfile
math::XYZVector XYZVector
int getmip()
Definition: HDShower.h:42
GflashProtonShowerProfile * theProtonProfile
const PreshowerLayer2Properties * layer2Properties(int onLayer2) const
Preshower Layer2 properties.
Definition: Calorimeter.cc:111
double getPathLengthAtShower()
Definition: GflashShowino.h:26
const std::map< CaloHitID, float > & getHits() override
not been done.
void setTrackParameters(const XYZNormal &normal, double X0depthoffset, const FSimTrack &theTrack)
void loadFromEcalEndcap(edm::PCaloHitContainer &c) const
void setID(const int id)
Definition: RawParticle.cc:102
const XYZTLorentzVector & deltaMom() const
Returns the actual energy lost.
std::vector< double > radiusPreshowerCorrections_
const HCALProperties * hcalProperties(int onHcal) const
HCAL properties.
Definition: Calorimeter.cc:87
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double getMIPfraction(double energy, double eta)
vec1 & getCorrHFhad()
Definition: HCALResponse.h:49
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:133
int onVFcal() const
Definition: FSimTrack.h:111
void setSpotEnergy(double e) override
Set the spot energy.
Definition: HcalHitMaker.h:28
std::vector< double > theTailIntervals_
std::vector< double > samplingHBHE_
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
void loadMuonSimTracks(edm::SimTrackContainer &m) const
unsigned int nTracks() const
Number of tracks.
Definition: FSimEvent.cc:36
JetCorrectorParameters corr
Definition: classes.h:5
HCALResponse * myHDResponse_
vec1 & getCorrHFem()
Definition: HCALResponse.h:48
CaloGeometryHelper * myCalorimeter_
std::vector< double > timeShiftHB_
int ietaAbs() const
get the absolute value of the cell ieta
Definition: HcalDetId.cc:98
const double intLength[kNumberCalorimeter]
std::vector< double > theCoreIntervals_
void getGflashTrajectoryPoint(GflashTrajectoryPoint &point, double s) const
edm::EventID id() const
Method to return the EventId.
Definition: FSimEvent.cc:27
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
Definition: DetId.h:18
void loadFromHcal(edm::PCaloHitContainer &c) const
GammaFunctionGenerator * aGammaGenerator
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
CLHEP::Hep3Vector Gflash3Vector
Definition: Gflash3Vector.h:6
std::vector< double > mipValues_
static std::vector< std::pair< int, float > > myZero_
bool null() const
is this a null id ?
Definition: DetId.h:45
std::vector< std::pair< CaloHitID, float > > EEMapping_
void setTkPosition(const math::XYZVectorD &pos)
Definition: SimTrack.h:40
double getPathLengthOnEcal()
Definition: GflashShowino.h:25
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.h:114
math::XYZVector XYZPoint
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
HLT enums.
std::vector< unsigned int > evtsToDebug_
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:86
bool preshowerPresent() const
const std::map< CaloHitID, float > & getHitsMap()
void setPulledPadSurvivalProbability(double val)
Definition: EcalHitMaker.h:132
const RawParticle & layer2Entrance() const
The particle at Preshower Layer 2.
Definition: FSimTrack.h:130
int onLayer1() const
Definition: FSimTrack.h:91
std::vector< double > k_e
std::vector< GflashHit > & getGflashHitList()
int onHcal() const
Definition: FSimTrack.h:106
const std::map< CaloHitID, float > & getHits() override
Definition: HcalHitMaker.h:37
void setTkMomentum(const math::XYZTLorentzVectorD &mom)
Definition: SimTrack.h:42
double eta() const
Definition: RawParticle.h:267
const RawParticle & hcalEntrance() const
The particle at HCAL entrance.
Definition: FSimTrack.h:136
std::vector< std::pair< CaloHitID, float > > HMapping_
std::vector< double > samplingHO_
void setPreshowerPresent(bool ps)
Definition: EcalHitMaker.h:137
FastHFShowerLibrary * theHFShowerLibrary
DetId getClosestCell(const XYZPoint &point, bool ecal, bool central) const
EnergyLossSimulator * energyLossSimulator() const
Return the Energy Loss engine.
std::unique_ptr< KKCorrectionFactors > ecalCorrection
HSParameters * myHSParameters_
std::vector< SimTrack > SimTrackContainer
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:288
bool setDepth(double, bool inCm=false)
set the depth in X0 or Lambda0 units depending on showerType
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
bool computeShower()
Definition: HDRShower.cc:47
math::XYZTLorentzVector XYZTLorentzVector
Definition: RawParticle.h:15
void print() const
print the FBaseSimEvent in an intelligible way
void setHcal(HcalHitMaker *const myHcal)
set the HCAL address
Definition: EMShower.cc:714
void setGrid(EcalHitMaker *const myGrid)
set the grid address
Definition: EMShower.h:61
void EMShowerSimulation(const FSimTrack &myTrack, RandomEngineAndDistribution const *)
std::vector< double > timeShiftHE_
FSimTrack & track(int id) const
Return track with given Id.
double getMaximumOfShower() const
get the depth of the centre of gravity of the shower(s)
Definition: EMShower.h:58
void loadFromEcalBarrel(edm::PCaloHitContainer &c) const