CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CalorimetryManager.cc
Go to the documentation of this file.
1 //Framework headers
4 
5 // Fast Simulation headers
20 //#include "FastSimulation/Utilities/interface/Histos.h"
30 // New headers for Muon Mip Simulation
33 // Muon Brem
35 
36 //Gflash Hadronic Model
44 // STL headers
45 #include <vector>
46 #include <iostream>
47 
48 //CMSSW headers
51 //#include "DataFormats/EcalDetId/interface/EcalDetId.h"
52 
54 
55 using namespace edm;
56 
59 
60 std::vector<std::pair<int,float> > CalorimetryManager::myZero_ = std::vector<std::pair<int,float> >
61 (1,std::pair<int,float>(0,0.));
62 
64  myCalorimeter_(0),
65  myHistos(0),
66  random(0),initialized_(false)
67 {;}
68 
70  const edm::ParameterSet& fastCalo,
71  const edm::ParameterSet& fastMuECAL,
72  const edm::ParameterSet& fastMuHCAL,
73  const edm::ParameterSet& parGflash,
74  const RandomEngine* engine)
75  :
76  mySimEvent(aSimEvent),
77  random(engine),initialized_(false),
78  theMuonEcalEffects(0), theMuonHcalEffects (0)
79 
80 {
81 
84 
85  //Gflash
86  theProfile = new GflashHadronShowerProfile(parGflash);
87  thePiKProfile = new GflashPiKShowerProfile(parGflash);
90 
91  readParameters(fastCalo);
92 
93 // EBMapping_.resize(62000,myZero_);
94 // EEMapping_.resize(20000,myZero_);
95 // HMapping_.resize(10000,myZero_);
96  EBMapping_.resize(62000);
97  EEMapping_.resize(20000);
98  HMapping_.resize(10000);
99  theDetIds_.resize(10000);
100 
101  unsigned s=(unfoldedMode_)?5:1;
102  for(unsigned ic=0;ic<62000;++ic)
103  {
104  EBMapping_[ic].reserve(s);
105  if(ic<20000)
106  EEMapping_[ic].reserve(s);
107  if(ic<10000)
108  HMapping_[ic].reserve(s);
109  }
110 
111 
112 
113  myHistos = 0;
114 #ifdef FAMOSDEBUG
116  myHistos->book("h10",140,-3.5,3.5,100,-0.5,99.5);
117  myHistos->book("h20",150,0,150.,100,-0.5,99.5);
118  myHistos->book("h100",140,-3.5,3.5,100,0,0.1);
119  myHistos->book("h110",140,-3.5,3.5,100,0,10.);
120  myHistos->book("h120",200,-5.,5.,100,0,0.5);
121 
122  myHistos->book("h200",300,0,3.,100,0.,35.);
123  myHistos->book("h210",720,-M_PI,M_PI,100,0,35.);
124  myHistos->book("h212",720,-M_PI,M_PI,100,0,35.);
125 
126  myHistos->bookByNumber("h30",0,7,300,-3.,3.,100,0.,35.);
127  myHistos->book("h310",75,-3.,3.,"");
128  myHistos->book("h400",100,-10.,10.,100,0.,35.);
129  myHistos->book("h410",720,-M_PI,M_PI);
130 #endif
131  myCalorimeter_ =
132  new CaloGeometryHelper(fastCalo);
133  myHDResponse_ =
134  new HCALResponse(fastCalo.getParameter<edm::ParameterSet>("HCALResponse"),
135  random);
136  myHSParameters_ =
137  new HSParameters(fastCalo.getParameter<edm::ParameterSet>("HSParameters"));
138 
139  // Material Effects for Muons in ECAL (only EnergyLoss implemented so far)
140 
141  if ( fastMuECAL.getParameter<bool>("PairProduction") ||
142  fastMuECAL.getParameter<bool>("Bremsstrahlung") ||
143  fastMuECAL.getParameter<bool>("MuonBremsstrahlung") ||
144  fastMuECAL.getParameter<bool>("EnergyLoss") ||
145  fastMuECAL.getParameter<bool>("MultipleScattering") )
146  theMuonEcalEffects = new MaterialEffects(fastMuECAL,random);
147 
148  // Material Effects for Muons in HCAL (only EnergyLoss implemented so far)
149 
150  if ( fastMuHCAL.getParameter<bool>("PairProduction") ||
151  fastMuHCAL.getParameter<bool>("Bremsstrahlung") ||
152  fastMuHCAL.getParameter<bool>("MuonBremsstrahlung") ||
153  fastMuHCAL.getParameter<bool>("EnergyLoss") ||
154  fastMuHCAL.getParameter<bool>("MultipleScattering") )
155  theMuonHcalEffects = new MaterialEffects(fastMuHCAL,random);
156 
157 }
158 
160 {
161  unsigned size=firedCellsEB_.size();
162  for(unsigned ic=0;ic<size;++ic)
163  {
164  EBMapping_[firedCellsEB_[ic]].clear();
165  }
166  firedCellsEB_.clear();
167 
168  size=firedCellsEE_.size();
169  for(unsigned ic=0;ic<size;++ic)
170  {
171  EEMapping_[firedCellsEE_[ic]].clear();
172  }
173  firedCellsEE_.clear();
174 
175  size=firedCellsHCAL_.size();
176  for(unsigned ic=0;ic<size;++ic)
177  {
178  HMapping_[firedCellsHCAL_[ic]].clear();
179  }
180  firedCellsHCAL_.clear();
181 
182  ESMapping_.clear();
183  muonSimTracks.clear();
184 }
185 
187 {
188 #ifdef FAMOSDEBUG
189  myHistos->put("Famos.root");
190 #endif
191  if(myCalorimeter_) delete myCalorimeter_;
192  if(myHDResponse_) delete myHDResponse_;
193 
196 
197  if ( theProfile ) delete theProfile;
198 }
199 
201 {
202  if(evtsToDebug_.size())
203  {
204  std::vector<unsigned int>::const_iterator itcheck=find(evtsToDebug_.begin(),evtsToDebug_.end(),mySimEvent->id().event());
205  debug_=(itcheck!=evtsToDebug_.end());
206  if(debug_)
207  mySimEvent->print();
208  }
209  // Clear the content of the calorimeters
210  if(!initialized_)
211  {
213  for(int subdetn=1;subdetn<=4;++subdetn)
214  {
215  const std::vector<DetId>& ids(geom->getValidDetIds(DetId::Hcal,subdetn));
216  for (std::vector<DetId>::const_iterator i=ids.begin(); i!=ids.end(); i++)
217  {
218  HcalDetId myDetId(*i);
219  unsigned hi=myDetId.hashed_index();
220  theDetIds_[hi]=myDetId;
221  }
222  }
223 
224  // Check if the preshower is really available
226  {
227  std::cout << " WARNING " << std::endl;
228  std::cout << " The preshower simulation has been turned on; but no preshower geometry is available " << std::endl;
229  std::cout << " Disabling the preshower simulation " << std::endl;
230  simulatePreshower_ = false;
231  }
232 
233  initialized_=true;
234  }
235  clean();
236 
237  LogInfo("FastCalorimetry") << "Reconstructing " << (int) mySimEvent->nTracks() << " tracks." << std::endl;
238  for( int fsimi=0; fsimi < (int) mySimEvent->nTracks() ; ++fsimi) {
239 
240  FSimTrack& myTrack = mySimEvent->track(fsimi);
241 
242  int pid = abs(myTrack.type());
243 
244  if (debug_) {
245  LogDebug("FastCalorimetry") << " ===> pid = " << pid << std::endl;
246  }
247 
248 
249  // Check that the particle hasn't decayed
250  if(myTrack.noEndVertex()) {
251  // Simulate energy smearing for photon and electrons
252  if ( pid == 11 || pid == 22 ) {
253 
254 
255  if ( myTrack.onEcal() )
256  EMShowerSimulation(myTrack);
257  else if ( myTrack.onVFcal() )
258  reconstructHCAL(myTrack);
259 
260  } // electron or photon
261  else if (pid==13)
262  {
263  MuonMipSimulation(myTrack);
264  }
265  // Simulate energy smearing for hadrons (i.e., everything
266  // but muons... and SUSY particles that deserve a special
267  // treatment.
268  else if ( pid < 1000000 ) {
269  if ( myTrack.onHcal() || myTrack.onVFcal() ) {
270  if(optionHDSim_ == 0 ) reconstructHCAL(myTrack);
271  else HDShowerSimulation(myTrack);
272  }
273  } // pid < 1000000
274  } // myTrack.noEndVertex()
275  } // particle loop
276  // LogInfo("FastCalorimetry") << " Number of hits (barrel)" << EBMapping_.size() << std::endl;
277  // LogInfo("FastCalorimetry") << " Number of hits (Hcal)" << HMapping_.size() << std::endl;
278  // std::cout << " Nombre de hit (endcap)" << EEMapping_.size() << std::endl;
279 
280 } // reconstruct
281 
282 // Simulation of electromagnetic showers in PS, ECAL, HCAL
284  std::vector<const RawParticle*> thePart;
285  double X0depth;
286  if (debug_) {
287  LogDebug("FastCalorimetry") << " EMShowerSimulation " <<myTrack << std::endl;
288  }
289 
290  // std::cout << " Simulating " << myTrack << std::endl;
291 
292  // The Particle at ECAL entrance
293  // std::cout << " Before ecalEntrance " << std::endl;
294  myPart = myTrack.ecalEntrance();
295 
296  // protection against infinite loop.
297  if ( myTrack.type() == 22 && myPart.e()<0.055) return;
298 
299 
300  // Barrel or Endcap ?
301  int onEcal = myTrack.onEcal();
302  int onHcal = myTrack.onHcal();
303  int onLayer1 = myTrack.onLayer1();
304  int onLayer2 = myTrack.onLayer2();
305 
306  // The entrance in ECAL
307  XYZPoint ecalentrance = myPart.vertex().Vect();
308 
309  // std::cout << " Ecal entrance " << ecalentrance << std::endl;
310 
311  // The preshower
312  PreshowerHitMaker * myPreshower = NULL ;
313  if(simulatePreshower_ && (onLayer1 || onLayer2))
314  {
315  XYZPoint layer1entrance,layer2entrance;
316  XYZVector dir1,dir2;
317  if(onLayer1)
318  {
319  layer1entrance = XYZPoint(myTrack.layer1Entrance().vertex().Vect());
320  dir1 = XYZVector(myTrack.layer1Entrance().Vect().Unit());
321  }
322  if(onLayer2)
323  {
324  layer2entrance = XYZPoint(myTrack.layer2Entrance().vertex().Vect());
325  dir2 = XYZVector(myTrack.layer2Entrance().Vect().Unit());
326  }
327  // std::cout << " Layer1entrance " << layer1entrance << std::endl;
328  // std::cout << " Layer2entrance " << layer2entrance << std::endl;
329  myPreshower = new PreshowerHitMaker(myCalorimeter_,
330  layer1entrance,
331  dir1,
332  layer2entrance,
333  dir2,
335  myPreshower->setMipEnergy(mipValues_[0],mipValues_[1]);
336  }
337 
338  // The ECAL Properties
340  showerparam(myCalorimeter_->ecalProperties(onEcal),
341  myCalorimeter_->hcalProperties(onHcal),
342  myCalorimeter_->layer1Properties(onLayer1),
343  myCalorimeter_->layer2Properties(onLayer2),
346  RCFactor_,
347  RTFactor_);
348 
349  // Photons : create an e+e- pair
350  if ( myTrack.type() == 22 ) {
351 
352  // Depth for the first e+e- pair creation (in X0)
353  X0depth = -log(random->flatShoot()) * (9./7.);
354 
355  // Initialization
356  double eMass = 0.000510998902;
357  double xe=0;
358  double xm=eMass/myPart.e();
359  double weight = 0.;
360 
361  // Generate electron energy between emass and eGamma-emass
362  do {
363  xe = random->flatShoot()*(1.-2.*xm) + xm;
364  weight = 1. - 4./3.*xe*(1.-xe);
365  } while ( weight < random->flatShoot() );
366 
367  // Protection agains infinite loop in Famos Shower
368  if ( myPart.e()*xe < 0.055 || myPart.e()*(1.-xe) < 0.055 ) {
369 
370  if ( myPart.e() > 0.055 ) thePart.push_back(&myPart);
371 
372  } else {
373 
374  myElec = (myPart) * xe;
375  myPosi = (myPart) * (1.-xe);
378  thePart.push_back(&myElec);
379  thePart.push_back(&myPosi);
380  }
381  // Electrons
382  } else {
383 
384  X0depth = 0.;
385  if ( myPart.e() > 0.055 ) thePart.push_back(&myPart);
386 
387  }
388 
389  // After the different protections, this shouldn't happen.
390  if(thePart.size()==0)
391  {
392  if(myPreshower==NULL) return;
393  delete myPreshower;
394  return;
395  }
396 
397  // find the most energetic particle
398  double maxEnergy=-1.;
399  for(unsigned ip=0;ip < thePart.size();++ip)
400  if(thePart[ip]->e() > maxEnergy) maxEnergy = thePart[ip]->e();
401 
402  // Initialize the Grid in ECAL
403  int size = gridSize_;
404  if(maxEnergy>100) size=11;
405 // if ( maxEnergy < threshold5x5 ) size = 5;
406 // if ( maxEnergy < threshold3x3 ) size = 3;
407 
408 
409  EMShower theShower(random,aGammaGenerator,&showerparam,&thePart);
410 
411  double maxShower = theShower.getMaximumOfShower();
412  if (maxShower > 20.) maxShower = 2.; // simple pivot-searching protection
413 
414  double depth((X0depth + maxShower) *
416  XYZPoint meanShower = ecalentrance + myPart.Vect().Unit()*depth;
417 
418  // if(onEcal!=1) return ;
419 
420  // The closest crystal
421  DetId pivot(myCalorimeter_->getClosestCell(meanShower, true, onEcal==1));
422 
423  if(pivot.subdetId() == 0) { // further protection against avbsence of pivot
424  edm::LogWarning("CalorimetryManager") << "Pivot for egamma e = " << myTrack.hcalEntrance().e() << " is not found at depth " << depth << " and meanShower coordinates = " << meanShower << std::endl;
425  return;
426  }
427 
428  EcalHitMaker myGrid(myCalorimeter_,ecalentrance,pivot,onEcal,size,0,random);
429  // ^^^^
430  // for EM showers
433 
434  //maximumdepth dependence of the radiusfactorbehindpreshower
435  //First tuning: Shilpi Jain (Mar-Apr 2010); changed after tuning - Feb-July - Shilpi Jain
436  /* **************
437  myGrid.setRadiusFactor(radiusFactor_);
438  if(onLayer1 || onLayer2)
439  {
440  float b = radiusPreshowerCorrections_[0];
441  float a = radiusFactor_*( 1.+radiusPreshowerCorrections_[1]*radiusPreshowerCorrections_[0] );
442  float maxdepth = X0depth+theShower.getMaximumOfShower();
443  float newRadiusFactor = radiusFactor_;
444  if(myPart.e()<=250.)
445  {
446  newRadiusFactor = a/(1.+b*maxdepth);
447  }
448  myGrid.setRadiusFactor(newRadiusFactor);
449  }
450  else // otherwise use the normal radius factor
451  {
452  myGrid.setRadiusFactor(radiusFactor_);
453  }
454  ************** */
455  if(myTrack.onEcal() == 2) // if on EE
456  {
457  if( (onLayer1 || onLayer2) && myPart.e()<=250.)
458  {
459  double maxdepth = X0depth+theShower.getMaximumOfShower();
460  double newRadiusFactor = radiusFactorEE_ * aTerm/(1.+bTerm*maxdepth);
461  myGrid.setRadiusFactor(newRadiusFactor);
462  }
463  else // otherwise use the normal radius factor
464  {
466  }
467  }//if(myTrack.onEcal() == 2)
468  else // else if on EB
469  {
471  }
472  //(end of) changed after tuning - Feb-July - Shilpi Jain
473 
475 
476  // The shower simulation
477  myGrid.setTrackParameters(myPart.Vect().Unit(),X0depth,myTrack);
478 
479 // std::cout << " PS ECAL GAP HCAL X0 " << myGrid.ps1TotalX0()+myGrid.ps2TotalX0() << " " << myGrid.ecalTotalX0();
480 // std::cout << " " << myGrid.ecalHcalGapTotalX0() << " " << myGrid.hcalTotalX0() << std::endl;
481 // std::cout << " PS ECAL GAP HCAL L0 " << myGrid.ps1TotalL0()+myGrid.ps2TotalL0() << " " << myGrid.ecalTotalL0();
482 // std::cout << " " << myGrid.ecalHcalGapTotalL0() << " " << myGrid.hcalTotalL0() << std::endl;
483 // std::cout << "ECAL-HCAL " << myTrack.momentum().eta() << " " << myGrid.ecalHcalGapTotalL0() << std::endl;
484 //
485 // std::cout << " Grid created " << std::endl;
486  if(myPreshower) theShower.setPreshower(myPreshower);
487 
488  HcalHitMaker myHcalHitMaker(myGrid,(unsigned)0);
489 
490  theShower.setGrid(&myGrid);
491  theShower.setHcal(&myHcalHitMaker);
492  theShower.compute();
493  //myHistos->fill("h502", myPart->eta(),myGrid.totalX0());
494 
495  // Save the hits !
496  std::map<uint32_t,float>::const_iterator mapitr;
497  std::map<uint32_t,float>::const_iterator endmapitr=myGrid.getHits().end();
498  for(mapitr=myGrid.getHits().begin();mapitr!=endmapitr;++mapitr)
499  {
500  if(onEcal==1)
501  {
502  updateMap(EBDetId(mapitr->first).hashedIndex(), mapitr->second,myTrack.id(),EBMapping_,firedCellsEB_);
503  }
504 
505  else if(onEcal==2)
506  updateMap(EEDetId(mapitr->first).hashedIndex(), mapitr->second,myTrack.id(),EEMapping_,firedCellsEE_);
507  // std::cout << " Adding " <<mapitr->first << " " << mapitr->second <<std::endl;
508  }
509 
510  // Now fill the HCAL hits
511  endmapitr=myHcalHitMaker.getHits().end();
512  for(mapitr=myHcalHitMaker.getHits().begin();mapitr!=endmapitr;++mapitr)
513  {
514  updateMap(HcalDetId(mapitr->first).hashed_index(),mapitr->second,myTrack.id(),HMapping_,firedCellsHCAL_);
515  // std::cout << " Adding " <<mapitr->first << " " << mapitr->second <<std::endl;
516  }
517 
518  // delete the preshower
519  if(myPreshower!=0)
520  {
521  endmapitr=myPreshower->getHits().end();
522  for(mapitr=myPreshower->getHits().begin();mapitr!=endmapitr;++mapitr)
523  {
524  updateMap(mapitr->first,mapitr->second,myTrack.id(),ESMapping_);
525  // std::cout << " Adding " <<mapitr->first << " " << mapitr->second <<std::endl;
526  }
527  delete myPreshower;
528  // std::cout << " Deleting myPreshower " << std::endl;
529  }
530 
531 }
532 
533 
534 
535 // Simulation of electromagnetic showers in VFCAL
537  if(debug_) {
538  XYZTLorentzVector moment = track.momentum();
539  std::cout << "FASTEnergyReconstructor::reconstructECAL - " << std::endl
540  << " eta " << moment.eta() << std::endl
541  << " phi " << moment.phi() << std::endl
542  << " et " << moment.Et() << std::endl;
543  }
544 
545  int hit;
546 
547  bool central=track.onEcal()==1;
548 
549  //Reconstruct only electrons and photons.
550 
551  //deal with different conventions
552  // ParticlePropagator 1 <-> Barrel
553  // 2 <-> EC
554  // whereas for Artur(this code):
555  // 0 <-> Barrel
556  // 1 <-> EC
557  // 2 <-> VF
558  XYZTLorentzVector trackPosition;
559  if( track.onEcal() ) {
560  hit=track.onEcal()-1;
561  trackPosition=track.ecalEntrance().vertex();
562  } else {
563  hit=2;
564  trackPosition=track.vfcalEntrance().vertex();
565  }
566 
567  double pathEta = trackPosition.eta();
568  double pathPhi = trackPosition.phi();
569  double EGen = track.ecalEntrance().e();
570 
571 
572  double e=0.;
573  double sigma=0.;
574  // if full simulation and in HF, but without showering anyway...
575  if(hit == 2 && optionHDSim_ == 2 ) {
576  std::pair<double,double> response =
577  myHDResponse_->responseHCAL(0, EGen, pathEta, 0); // last par.= 0 = e/gamma
578  e = response.first;
579  sigma = response.second;
580  }
581 
582  double emeas = 0.;
583 
584  if(sigma>0.)
585  emeas = random->gaussShoot(e,sigma);
586 
587 
588  if(debug_)
589  std::cout << "FASTEnergyReconstructor::reconstructECAL : "
590  << " on-calo eta, phi = " << pathEta << " " << pathPhi << std::endl
591  << " Egen = " << EGen << std::endl
592  << " Eres = " << e << std::endl
593  << " sigma = " << sigma << std::endl
594  << " Emeas = " << emeas << std::endl;
595 
596 
597  if(debug_)
598  std::cout << "FASTEnergyReconstructor::reconstructECAL : "
599  << " Track position - " << trackPosition.Vect()
600  << " bool central - " << central
601  << " hit - " << hit << std::endl;
602 
603  DetId detid;
604  if( hit==2 )
605  detid = myCalorimeter_->getClosestCell(trackPosition.Vect(),false,central);
606  // Check that the detid is HCAL forward
607  HcalDetId hdetid(detid);
608  if(!hdetid.subdetId()!=HcalForward) return;
609 
610  if(debug_)
611  std::cout << "FASTEnergyReconstructor::reconstructECAL : "
612  << " CellID - " << detid.rawId() << std::endl;
613 
614  if( hit != 2 || emeas > 0.)
615  if(!detid.null())
616  {
617  updateMap(hdetid.hashed_index(),emeas,track.id(),HMapping_,firedCellsHCAL_);
618  }
619 
620 }
621 
622 
624 {
625  int hit;
626  int pid = abs(myTrack.type());
627  if (debug_) {
628  LogDebug("FastCalorimetry") << " reconstructHCAL " << myTrack << std::endl;
629  }
630 
631  // FSimTrack myTrack = mySimEvent.track(fsimi);
632 
633  // int pid=abs(myTrack.type());
634  // std::cout << "reconstructHCAL " << std::endl;
635 
636  XYZTLorentzVector trackPosition;
637  if (myTrack.onHcal()) {
638  trackPosition=myTrack.hcalEntrance().vertex();
639  hit = myTrack.onHcal()-1;
640  } else {
641  trackPosition=myTrack.vfcalEntrance().vertex();
642  hit = 2;
643  }
644 
645  double pathEta = trackPosition.eta();
646  double pathPhi = trackPosition.phi();
647  // double pathTheta = trackPosition.theta();
648 
649  double EGen = myTrack.hcalEntrance().e();
650  double e = 0.;
651  double sigma = 0.;
652  double emeas = 0.;
653 
654  if(pid == 13) {
655  // std::cout << " We should not be here " << std::endl;
656  std::pair<double,double> response =
657  myHDResponse_->responseHCAL(0, EGen, pathEta, 2); // 2=muon
658  emeas = response.first;
659  if(debug_)
660  LogDebug("FastCalorimetry") << "CalorimetryManager::reconstructHCAL - MUON !!!" << std::endl;
661  }
662  else if( pid == 22 || pid == 11)
663  {
664 
665  std::pair<double,double> response =
666  myHDResponse_->responseHCAL(0, EGen, pathEta, 0); // last par. = 0 = e/gamma
667  e = response.first; //
668  sigma = response.second; //
669  emeas = random->gaussShoot(e,sigma); //
670 
671  // cout << "CalorimetryManager::reconstructHCAL - e/gamma !!!" << std::endl;
672  if(debug_)
673  LogDebug("FastCalorimetry") << "CalorimetryManager::reconstructHCAL - e/gamma !!!" << std::endl;
674  }
675  else {
676  e = myHDResponse_->getHCALEnergyResponse(EGen,hit);
677  sigma = myHDResponse_->getHCALEnergyResolution(EGen, hit);
678 
679  emeas = random->gaussShoot(e,sigma);
680  }
681 
682 
683  if(debug_)
684  LogDebug("FastCalorimetry") << "CalorimetryManager::reconstructHCAL - on-calo "
685  << " eta = " << pathEta
686  << " phi = " << pathPhi
687  << " Egen = " << EGen
688  << " Eres = " << e
689  << " sigma = " << sigma
690  << " Emeas = " << emeas << std::endl;
691 
692  if(emeas > 0.) {
693  DetId cell = myCalorimeter_->getClosestCell(trackPosition.Vect(),false,false);
694  updateMap(HcalDetId(cell).hashed_index(), emeas, myTrack.id(),HMapping_,firedCellsHCAL_);
695  }
696 }
697 
699 {
700  // TimeMe t(" FASTEnergyReconstructor::HDShower");
701  XYZTLorentzVector moment = myTrack.momentum();
702 
703  if(debug_)
704  LogDebug("FastCalorimetry")
705  << "CalorimetryManager::HDShowerSimulation - track param."
706  << std::endl
707  << " eta = " << moment.eta() << std::endl
708  << " phi = " << moment.phi() << std::endl
709  << " et = " << moment.Et() << std::endl
710  << " e = " << myTrack.hcalEntrance().e() << std::endl;
711 
712  if (debug_) {
713  LogDebug("FastCalorimetry") << " HDShowerSimulation " << myTrack << std::endl;
714  }
715 
716 
717  int hit;
718  // int pid = abs(myTrack.type());
719 
720  XYZTLorentzVector trackPosition;
721  if ( myTrack.onEcal() ) {
722  trackPosition=myTrack.ecalEntrance().vertex();
723  hit = myTrack.onEcal()-1; //
724  myPart = myTrack.ecalEntrance();
725  } else if ( myTrack.onVFcal()) {
726  trackPosition=myTrack.vfcalEntrance().vertex();
727  hit = 2;
728  myPart = myTrack.vfcalEntrance();
729  }
730  else
731  {
732  LogDebug("FastCalorimetry") << " The particle is not in the acceptance " << std::endl;
733  return;
734  }
735 
736  // int onHCAL = hit + 1; - specially for myCalorimeter->hcalProperties(onHCAL)
737  // (below) to get VFcal properties ...
738  int onHCAL = hit + 1;
739  int onECAL = myTrack.onEcal();
740 
741  double pathEta = trackPosition.eta();
742  double pathPhi = trackPosition.phi();
743  // double pathTheta = trackPosition.theta();
744 
745  double eint = moment.e();
746  double eGen = myTrack.hcalEntrance().e();
747  double e = 0.;
748  double sigma = 0.;
749 
750  double emeas = 0.;
751 
752  //===========================================================================
753  if(eGen > 0.) {
754 
755  // ECAL and HCAL properties to get
757  theHDShowerparam(myCalorimeter_->ecalProperties(onECAL),
760 
761  //Making ECAL Grid (and segments calculation)
762  XYZPoint caloentrance;
763  XYZVector direction;
764  if(myTrack.onEcal())
765  {
766  caloentrance = myTrack.ecalEntrance().vertex().Vect();
767  direction = myTrack.ecalEntrance().Vect().Unit();
768  }
769  else if(myTrack.onHcal())
770  {
771  caloentrance = myTrack.hcalEntrance().vertex().Vect();
772  direction = myTrack.hcalEntrance().Vect().Unit();
773  }
774  else
775  {
776  caloentrance = myTrack.vfcalEntrance().vertex().Vect();
777  direction = myTrack.vfcalEntrance().Vect().Unit();
778  }
779 
780  if(debug_)
781  LogDebug("FastCalorimetry")
782  << "CalorimetryManager::HDShowerSimulation - on-calo 1 "
783  << std::endl
784  << " onEcal = " << myTrack.onEcal() << std::endl
785  << " onHcal = " << myTrack.onHcal() << std::endl
786  << " onVFcal = " << myTrack.onVFcal() << std::endl
787  << " position = " << caloentrance << std::endl;
788 
789 
790  DetId pivot;
791  if(myTrack.onEcal())
792  {
793  pivot=myCalorimeter_->getClosestCell(caloentrance,
794  true, myTrack.onEcal()==1);
795  }
796  else if(myTrack.onHcal())
797  {
798  // std::cout << " CalorimetryManager onHcal " << myTrack.onHcal() << " caloentrance" << caloentrance << std::endl;
799  pivot=myCalorimeter_->getClosestCell(caloentrance,
800  false, false);
801  }
802 
803  EcalHitMaker myGrid(myCalorimeter_,caloentrance,pivot,
804  pivot.null()? 0 : myTrack.onEcal(),hdGridSize_,1,
805  random);
806  // 1=HAD shower
807 
808  myGrid.setTrackParameters(direction,0,myTrack);
809  // Build the FAMOS HCAL
810  HcalHitMaker myHcalHitMaker(myGrid,(unsigned)1);
811 
812  // Shower simulation
813  bool status = false;
814  int mip = 2;
815  // Use HFShower for HF
816  if ( !myTrack.onEcal() && !myTrack.onHcal() ) {
817  // std::cout << "CalorimetryManager::HDShowerSimulation(): track entrance = "
818  // << myTrack.vfcalEntrance().vertex().X() << " "
819  // << myTrack.vfcalEntrance().vertex().Y() << " "
820  // << myTrack.vfcalEntrance().vertex().Z() << " "
821  // << " , Energy (Gen/Scale) = " << eGen << " " << e << std::endl;
822 
823  // Warning : We give here the particle energy with the response
824  // but without the resolution/gaussian smearing
825  // For HF, the resolution is due to the PE statistic
826 
827  HFShower theShower(random,
828  &theHDShowerparam,
829  &myGrid,
830  &myHcalHitMaker,
831  onECAL,
832  eGen);
833  // eGen);
834  // e); // PV Warning : temporarly set the energy to the generated E
835 
836  status = theShower.compute();
837  } else {
838  if(hdSimMethod_ == 0) {
839  HDShower theShower(random,
840  &theHDShowerparam,
841  &myGrid,
842  &myHcalHitMaker,
843  onECAL,
844  eGen);
845  status = theShower.compute();
846  mip = theShower.getmip();
847  }
848  else if (hdSimMethod_ == 1) {
849  HDRShower theShower(random,
850  &theHDShowerparam,
851  &myGrid,
852  &myHcalHitMaker,
853  onECAL,
854  eGen);
855  status = theShower.computeShower();
856  mip = 2;
857  }
858  else if (hdSimMethod_ == 2 ) {
859  // std::cout << "Using GflashHadronShowerProfile hdSimMethod_ == 2" << std::endl;
860 
861  //dynamically loading a corresponding profile by the particle type
862  int particleType = myTrack.type();
864  if(particleType == -2212) theProfile = theAntiProtonProfile;
865  else if(particleType == 2212) theProfile = theProtonProfile;
866 
867  //input variables for GflashHadronShowerProfile
868  int showerType = 99 + myTrack.onEcal();
869  double globalTime = 150.0; // a temporary reference hit time in nanosecond
870  float charge = (float)(myTrack.charge());
871  Gflash3Vector gfpos(trackPosition.X(),trackPosition.Y(),trackPosition.Z());
872  Gflash3Vector gfmom(moment.X(),moment.Y(),moment.Z());
873 
874  theProfile->initialize(showerType,eGen,globalTime,charge,gfpos,gfmom);
877 
878  //make hits
879  std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
880  std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
881  std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
882 
884  bool result;
885 
886  for( ; spotIter != spotIterEnd; spotIter++){
887 
888  double pathLength = theProfile->getGflashShowino()->getPathLengthAtShower()
889  + (30*100/eGen)*(spotIter->getTime() - globalTime);
890 
891  double currentDepth = std::max(0.0,pathLength - theProfile->getGflashShowino()->getPathLengthOnEcal());
892 
893  //find the the showino position at the currentDepth
894  GflashTrajectoryPoint trajectoryPoint;
895  theProfile->getGflashShowino()->getHelix()->getGflashTrajectoryPoint(trajectoryPoint,pathLength);
896  Gflash3Vector positionAtCurrentDepth = trajectoryPoint.getPosition();
897  //find radial distrance
898  Gflash3Vector lateralDisplacement = positionAtCurrentDepth - spotIter->getPosition()/CLHEP::cm;
899  double rShower = lateralDisplacement.r();
900  double azimuthalAngle = lateralDisplacement.phi();
901 
902  whichCalor = Gflash::getCalorimeterNumber(positionAtCurrentDepth);
903 
904  if(whichCalor==Gflash::kESPM || whichCalor==Gflash::kENCA) {
905  bool statusPad = myGrid.getPads(currentDepth,true);
906  if(!statusPad) continue;
907  myGrid.setSpotEnergy(1.2*spotIter->getEnergy()/CLHEP::GeV);
908  result = myGrid.addHit(rShower/Gflash::intLength[Gflash::kESPM],azimuthalAngle,0);
909  }
910  else if(whichCalor==Gflash::kHB || whichCalor==Gflash::kHE) {
911  bool setHDdepth = myHcalHitMaker.setDepth(currentDepth,true);
912  if(!setHDdepth) continue;
913  myHcalHitMaker.setSpotEnergy(1.4*spotIter->getEnergy()/CLHEP::GeV);
914  result = myHcalHitMaker.addHit(rShower/Gflash::intLength[Gflash::kHB],azimuthalAngle,0);
915  }
916  }
917  status = true;
918  }
919  else {
920  edm::LogInfo("FastSimulationCalorimetry") << " SimMethod " << hdSimMethod_ <<" is NOT available ";
921  }
922  }
923 
924 
925  if(status) {
926 
927  // Here to switch between simple formulae and parameterized response
928  if(optionHDSim_ == 1) {
929  e = myHDResponse_->getHCALEnergyResponse (eGen, hit);
930  sigma = myHDResponse_->getHCALEnergyResolution(eGen, hit);
931  }
932  else { // optionHDsim == 2
933  std::pair<double,double> response =
934  myHDResponse_->responseHCAL(mip, eGen, pathEta, 1); // 1=hadron
935  e = response.first;
936  sigma = response.second;
937  }
938 
939  emeas = random->gaussShoot(e,sigma);
940  double correction = emeas / eGen;
941 
942  // RespCorrP factors (ECAL and HCAL separately) calculation
943  respCorr(eint);
944 
945  if(debug_)
946  LogDebug("FastCalorimetry")
947  << "CalorimetryManager::HDShowerSimulation - on-calo 2" << std::endl
948  << " eta = " << pathEta << std::endl
949  << " phi = " << pathPhi << std::endl
950  << " Egen = " << eGen << std::endl
951  << " Eres = " << e << std::endl
952  << " sigma = " << sigma << std::endl
953  << " Emeas = " << emeas << std::endl
954  << " corr = " << correction << std::endl
955  << " mip = " << mip << std::endl;
956 
957 
958  // was map<unsigned,double> but CaloHitMaker uses float
959  std::map<unsigned,float>::const_iterator mapitr;
960  std::map<unsigned,float>::const_iterator endmapitr;
961  if(myTrack.onEcal() > 0) {
962  // Save ECAL hits
963  endmapitr=myGrid.getHits().end();
964  for(mapitr=myGrid.getHits().begin(); mapitr!=endmapitr; ++mapitr) {
965  double energy = mapitr->second;
966  energy *= correction; // RESCALING
967  energy *= ecorr;
968 
969  if(energy > 0.000001) {
970  if(onECAL==1)
971  updateMap(EBDetId(mapitr->first).hashedIndex(),energy,myTrack.id(),EBMapping_,firedCellsEB_);
972 
973  else if(onECAL==2)
974  updateMap(EEDetId(mapitr->first).hashedIndex(),energy,myTrack.id(),EEMapping_,firedCellsEE_);
975 
976  if(debug_)
977  LogDebug("FastCalorimetry") << " ECAL cell " << mapitr->first << " added, E = "
978  << energy << std::endl;
979  }
980  }
981  }
982 
983  // Save HCAL hits
984  endmapitr=myHcalHitMaker.getHits().end();
985  for(mapitr=myHcalHitMaker.getHits().begin(); mapitr!=endmapitr; ++mapitr) {
986  double energy = mapitr->second;
987  energy *= correction; // RESCALING
988  energy *= hcorr;
989 
990  updateMap(HcalDetId(mapitr->first).hashed_index(),energy,myTrack.id(),HMapping_,firedCellsHCAL_);
991  if(debug_)
992  LogDebug("FastCalorimetry") << " HCAL cell "
993  << mapitr->first << " added E = "
994  << mapitr->second << std::endl;
995  }
996  }
997  else { // shower simulation failed
998 // std::cout << " Shower simulation failed " << trackPosition.Vect() << std::endl;
999 // std::cout << " The FSimTrack " << myTrack << std::endl;
1000 // std::cout << " HF entrance on VFcal" << myTrack.onVFcal() << std::endl;
1001 // std::cout << " trackPosition.eta() " << trackPosition.eta() << std::endl;
1002  if(myTrack.onHcal() || myTrack.onVFcal())
1003  {
1004  DetId cell = myCalorimeter_->getClosestCell(trackPosition.Vect(),false,false);
1005  updateMap(HcalDetId(cell).hashed_index(),emeas,myTrack.id(),HMapping_,firedCellsHCAL_);
1006  if(debug_)
1007  LogDebug("FastCalorimetry") << " HCAL simple cell "
1008  << cell.rawId() << " added E = "
1009  << emeas << std::endl;
1010  }
1011  }
1012 
1013  } // e > 0. ...
1014 
1015  if(debug_)
1016  LogDebug("FastCalorimetry") << std::endl << " FASTEnergyReconstructor::HDShowerSimulation finished "
1017  << std::endl;
1018 }
1019 
1020 
1022 {
1023  // TimeMe t(" FASTEnergyReconstructor::HDShower");
1024  XYZTLorentzVector moment = myTrack.momentum();
1025 
1026  // Backward compatibility behaviour
1027  if(!theMuonHcalEffects)
1028  {
1029  if(myTrack.onHcal() || myTrack.onVFcal() )
1030  reconstructHCAL(myTrack);
1031 
1032  return;
1033  }
1034 
1035  if(debug_)
1036  LogDebug("FastCalorimetry") << "CalorimetryManager::MuonMipSimulation - track param."
1037  << std::endl
1038  << " eta = " << moment.eta() << std::endl
1039  << " phi = " << moment.phi() << std::endl
1040  << " et = " << moment.Et() << std::endl;
1041 
1042  int hit;
1043  // int pid = abs(myTrack.type());
1044 
1045  XYZTLorentzVector trackPosition;
1046  if ( myTrack.onEcal() ) {
1047  trackPosition=myTrack.ecalEntrance().vertex();
1048  hit = myTrack.onEcal()-1; //
1049  myPart = myTrack.ecalEntrance();
1050  } else if ( myTrack.onVFcal()) {
1051  trackPosition=myTrack.vfcalEntrance().vertex();
1052  hit = 2;
1053  myPart = myTrack.vfcalEntrance();
1054  }
1055  else
1056  {
1057  LogDebug("FastCalorimetry") << " The particle is not in the acceptance " << std::endl;
1058  return;
1059  }
1060 
1061  // int onHCAL = hit + 1; - specially for myCalorimeter->hcalProperties(onHCAL)
1062  // (below) to get VFcal properties ...
1063  // not needed ?
1064  // int onHCAL = hit + 1;
1065  int onECAL = myTrack.onEcal();
1066 
1067  // double pathEta = trackPosition.eta();
1068  // double pathPhi = trackPosition.phi();
1069  // double pathTheta = trackPosition.theta();
1070 
1071  //===========================================================================
1072 
1073  // ECAL and HCAL properties to get
1074 
1075  //Making ECAL Grid (and segments calculation)
1076  XYZPoint caloentrance;
1077  XYZVector direction;
1078  if(myTrack.onEcal())
1079  {
1080  caloentrance = myTrack.ecalEntrance().vertex().Vect();
1081  direction = myTrack.ecalEntrance().Vect().Unit();
1082  }
1083  else if(myTrack.onHcal())
1084  {
1085  caloentrance = myTrack.hcalEntrance().vertex().Vect();
1086  direction = myTrack.hcalEntrance().Vect().Unit();
1087  }
1088  else
1089  {
1090  caloentrance = myTrack.vfcalEntrance().vertex().Vect();
1091  direction = myTrack.vfcalEntrance().Vect().Unit();
1092  }
1093 
1094  DetId pivot;
1095  if(myTrack.onEcal())
1096  {
1097  pivot=myCalorimeter_->getClosestCell(caloentrance,
1098  true, myTrack.onEcal()==1);
1099  }
1100  else if(myTrack.onHcal())
1101  {
1102  // std::cout << " CalorimetryManager onHcal " << myTrack.onHcal() << " caloentrance" << caloentrance << std::endl;
1103  pivot=myCalorimeter_->getClosestCell(caloentrance,
1104  false, false);
1105  }
1106 
1107  EcalHitMaker myGrid(myCalorimeter_,caloentrance,pivot,
1108  pivot.null()? 0 : myTrack.onEcal(),hdGridSize_,0,
1109  random);
1110  // 0 =EM shower -> Unit = X0
1111 
1112  myGrid.setTrackParameters(direction,0,myTrack);
1113 
1114  // Now get the path in the Preshower, ECAL and HCAL along a straight line extrapolation
1115  // but only those in the ECAL are used
1116 
1117  const std::vector<CaloSegment>& segments=myGrid.getSegments();
1118  unsigned nsegments=segments.size();
1119 
1120  int ifirstHcal=-1;
1121  int ilastEcal=-1;
1122  EnergyLossSimulator* energyLossECAL = 0;
1123  //New muon brem effect in Calorimetry 16th Jan 2011 (Sandro Fonseca de Souza and Andre Sznajder UERJ/Brazil)
1124  MuonBremsstrahlungSimulator* muonBremECAL = 0;
1125 
1127  // Muon brem in ECAL
1129 
1130  for(unsigned iseg=0;iseg<nsegments&&ifirstHcal<0;++iseg)
1131  {
1132 
1133  // in the ECAL, there are two types of segments: PbWO4 and GAP
1134  float segmentSizeinX0=segments[iseg].X0length();
1135 
1136  // Martijn - insert your computations here
1137  float energy=0.0;
1138  if (segmentSizeinX0>0.001 && segments[iseg].material()==CaloSegment::PbWO4 ) {
1139  // The energy loss simulator
1140  float charge = (float)(myTrack.charge());
1141  ParticlePropagator theMuon(moment,trackPosition,charge,0);
1142  theMuon.setID(-(int)charge*13);
1143  if ( energyLossECAL ) {
1144  energyLossECAL->updateState(theMuon, segmentSizeinX0);
1145  energy = energyLossECAL->deltaMom().E();
1146  moment -= energyLossECAL->deltaMom();
1147  }
1148  }
1149  // that's all for ECAL, Florian
1150  // Save the hit only if it is a crystal
1151  if(segments[iseg].material()==CaloSegment::PbWO4)
1152  {
1153  myGrid.getPads(segments[iseg].sX0Entrance()+segmentSizeinX0*0.5);
1154  myGrid.setSpotEnergy(energy);
1155  myGrid.addHit(0.,0.);
1156  ilastEcal=iseg;
1157  }
1158  // Check for end of loop:
1159  if(segments[iseg].material()==CaloSegment::HCAL)
1160  {
1161  ifirstHcal=iseg;
1162  }
1163  }
1164 
1165  // Position of Ecal Exit
1166  math::XYZVector ecalExit;
1167  if(ilastEcal>=0)
1168  math::XYZVector ecalExit=segments[ilastEcal].exit();
1169 
1170  // Position of HCAL entrance
1171  math::XYZVector hcalEntrance;
1172  if(ifirstHcal>=0)
1173  hcalEntrance=segments[ifirstHcal].entrance();
1174 
1175  // dummy HCAL exit, just to test the FSimTrack saving mechanism
1176  math::XYZVector hcalExit;
1177  if(ifirstHcal>=0)
1178  hcalExit=segments[ifirstHcal].exit();
1179 
1180 
1181  // Build the FAMOS HCAL
1182  HcalHitMaker myHcalHitMaker(myGrid,(unsigned)2);
1183  // float mipenergy=0.1;
1184  // Create the helix with the stepping helix propagator
1185  // to add a hit, just do
1186  // myHcalHitMaker.setSpotEnergy(mipenergy);
1187  // myHcalHitMaker.addHit(hcalEntrance);
1191  int ilastHcal=-1;
1192  float mipenergy=0.0;
1193  EnergyLossSimulator* energyLossHCAL = 0;
1194  // Implementation of Muon Brem Effect in HCAL (Sandro Fonseca de Souza and Andre Sznajder UERJ/Brazil)
1195  //Date 01/16/2011
1196  MuonBremsstrahlungSimulator* muonBremHCAL = 0;
1197 
1199  // Muon Brem effect
1201 
1202  if(ifirstHcal>0 && energyLossHCAL){
1203  for(unsigned iseg=ifirstHcal;iseg<nsegments;++iseg)
1204  {
1205  float segmentSizeinX0=segments[iseg].X0length();
1206  if (segmentSizeinX0>0.001 && segments[iseg].material()==CaloSegment::HCAL ) {
1207  // The energy loss simulator
1208  float charge = (float)(myTrack.charge());
1209  ParticlePropagator theMuon(moment,trackPosition,charge,0);
1210  theMuon.setID(-(int)charge*13);
1211  energyLossHCAL->updateState(theMuon, segmentSizeinX0);
1212  mipenergy = energyLossHCAL->deltaMom().E();
1213  moment -= energyLossHCAL->deltaMom();
1214  myHcalHitMaker.setSpotEnergy(mipenergy);
1215  myHcalHitMaker.addHit(segments[iseg].entrance());
1216  }
1217  if(segments[iseg].material()==CaloSegment::HCAL)
1218  {
1219  ilastHcal=iseg;
1220  }
1221  }
1222  }
1223  if(ilastHcal>=0) // Update hcalExit position from 'dummy' HCAL entrance to 'staight line' HCAL exit
1224  hcalExit=segments[ilastHcal].exit();
1225 
1230  //
1231 
1232 
1233 
1234  // Copy the muon SimTrack (Only for Energy loss)
1235  FSimTrack muonTrack(myTrack);
1236  if(energyLossHCAL) {
1237  muonTrack.setTkPosition(hcalExit);
1238  muonTrack.setTkMomentum(moment);
1239  } else if(energyLossECAL) {
1240  muonTrack.setTkPosition(ecalExit);
1241  muonTrack.setTkMomentum(moment);
1242  } // else just leave tracker surface position and momentum...
1243 
1244  muonSimTracks.push_back(muonTrack);
1245 
1246 
1247  // no need to change below this line
1248  std::map<unsigned,float>::const_iterator mapitr;
1249  std::map<unsigned,float>::const_iterator endmapitr;
1250  if(myTrack.onEcal() > 0) {
1251  // Save ECAL hits
1252  endmapitr=myGrid.getHits().end();
1253  for(mapitr=myGrid.getHits().begin(); mapitr!=endmapitr; ++mapitr) {
1254  double energy = mapitr->second;
1255  if(onECAL==1)
1256  {
1257  updateMap(EBDetId(mapitr->first).hashedIndex(),energy,myTrack.id(),EBMapping_,firedCellsEB_);
1258  }
1259  else if(onECAL==2)
1260  {
1261  updateMap(EEDetId(mapitr->first).hashedIndex(),energy,myTrack.id(),EEMapping_,firedCellsEE_);
1262  }
1263 
1264  if(debug_)
1265  LogDebug("FastCalorimetry") << " ECAL cell " << mapitr->first << " added, E = "
1266  << energy << std::endl;
1267  }
1268  }
1269 
1270  // Save HCAL hits
1271  endmapitr=myHcalHitMaker.getHits().end();
1272  for(mapitr=myHcalHitMaker.getHits().begin(); mapitr!=endmapitr; ++mapitr) {
1273  double energy = mapitr->second;
1274  {
1275  updateMap(HcalDetId(mapitr->first).hashed_index(),energy,myTrack.id(),HMapping_,firedCellsHCAL_);
1276  }
1277  if(debug_)
1278  LogDebug("FastCalorimetry") << " HCAL cell "
1279  << mapitr->first << " added E = "
1280  << mapitr->second << std::endl;
1281  }
1282 
1283  if(debug_)
1284  LogDebug("FastCalorimetry") << std::endl << " FASTEnergyReconstructor::MipShowerSimulation finished "
1285  << std::endl;
1286 }
1287 
1288 
1290 
1291  edm::ParameterSet ECALparameters = fastCalo.getParameter<edm::ParameterSet>("ECAL");
1292 
1293  evtsToDebug_ = fastCalo.getUntrackedParameter<std::vector<unsigned int> >("EvtsToDebug",std::vector<unsigned>());
1294  debug_ = fastCalo.getUntrackedParameter<bool>("Debug");
1295 
1296  gridSize_ = ECALparameters.getParameter<int>("GridSize");
1297  spotFraction_ = ECALparameters.getParameter<double>("SpotFraction");
1298  pulledPadSurvivalProbability_ = ECALparameters.getParameter<double>("FrontLeakageProbability");
1299  crackPadSurvivalProbability_ = ECALparameters.getParameter<double>("GapLossProbability");
1300  theCoreIntervals_ = ECALparameters.getParameter<std::vector<double> >("CoreIntervals");
1301  theTailIntervals_ = ECALparameters.getParameter<std::vector<double> >("TailIntervals");
1302 
1303  RCFactor_ = ECALparameters.getParameter<double>("RCFactor");
1304  RTFactor_ = ECALparameters.getParameter<double>("RTFactor");
1305  //changed after tuning - Feb-July - Shilpi Jain
1306  // radiusFactor_ = ECALparameters.getParameter<double>("RadiusFactor");
1307  radiusFactorEE_ = ECALparameters.getParameter<double>("RadiusFactorEE");
1308  radiusFactorEB_ = ECALparameters.getParameter<double>("RadiusFactorEB");
1309  //(end of) changed after tuning - Feb-July - Shilpi Jain
1310  radiusPreshowerCorrections_ = ECALparameters.getParameter<std::vector<double> >("RadiusPreshowerCorrections");
1311  aTerm = 1.+radiusPreshowerCorrections_[1]*radiusPreshowerCorrections_[0];
1312  bTerm = radiusPreshowerCorrections_[0];
1313  mipValues_ = ECALparameters.getParameter<std::vector<double> >("MipsinGeV");
1314  simulatePreshower_ = ECALparameters.getParameter<bool>("SimulatePreshower");
1315 
1316  if(gridSize_ <1) gridSize_= 7;
1317  if(pulledPadSurvivalProbability_ <0. || pulledPadSurvivalProbability_>1 ) pulledPadSurvivalProbability_= 1.;
1318  if(crackPadSurvivalProbability_ <0. || crackPadSurvivalProbability_>1 ) crackPadSurvivalProbability_= 0.9;
1319 
1320  LogInfo("FastCalorimetry") << " Fast ECAL simulation parameters " << std::endl;
1321  LogInfo("FastCalorimetry") << " =============================== " << std::endl;
1322  if(simulatePreshower_)
1323  LogInfo("FastCalorimetry") << " The preshower is present " << std::endl;
1324  else
1325  LogInfo("FastCalorimetry") << " The preshower is NOT present " << std::endl;
1326  LogInfo("FastCalorimetry") << " Grid Size : " << gridSize_ << std::endl;
1327  if(spotFraction_>0.)
1328  LogInfo("FastCalorimetry") << " Spot Fraction : " << spotFraction_ << std::endl;
1329  else
1330  {
1331  LogInfo("FastCalorimetry") << " Core of the shower " << std::endl;
1332  for(unsigned ir=0; ir < theCoreIntervals_.size()/2;++ir)
1333  {
1334  LogInfo("FastCalorimetry") << " r < " << theCoreIntervals_[ir*2] << " R_M : " << theCoreIntervals_[ir*2+1] << " ";
1335  }
1336  LogInfo("FastCalorimetry") << std::endl;
1337 
1338  LogInfo("FastCalorimetry") << " Tail of the shower " << std::endl;
1339  for(unsigned ir=0; ir < theTailIntervals_.size()/2;++ir)
1340  {
1341  LogInfo("FastCalorimetry") << " r < " << theTailIntervals_[ir*2] << " R_M : " << theTailIntervals_[ir*2+1] << " ";
1342  }
1343  //changed after tuning - Feb-July - Shilpi Jain
1344  // LogInfo("FastCalorimetry") << "Radius correction factor " << radiusFactor_ << std::endl;
1345  LogInfo("FastCalorimetry") << "Radius correction factors: EB & EE " << radiusFactorEB_ << " : "<< radiusFactorEE_ << std::endl;
1346  //(end of) changed after tuning - Feb-July - Shilpi Jain
1347  LogInfo("FastCalorimetry") << std::endl;
1348  if(mipValues_.size()>2) {
1349  LogInfo("FastCalorimetry") << "Improper number of parameters for the preshower ; using 95keV" << std::endl;
1350  mipValues_.clear();
1351  mipValues_.resize(2,0.000095);
1352  }
1353  }
1354 
1355  LogInfo("FastCalorimetry") << " FrontLeakageProbability : " << pulledPadSurvivalProbability_ << std::endl;
1356  LogInfo("FastCalorimetry") << " GapLossProbability : " << crackPadSurvivalProbability_ << std::endl;
1357 
1358 
1359  // RespCorrP: p (momentum), ECAL and HCAL corrections = f(p)
1360  edm::ParameterSet CalorimeterParam = fastCalo.getParameter<edm::ParameterSet>("CalorimeterProperties");
1361 
1362  rsp = CalorimeterParam.getParameter<std::vector<double> >("RespCorrP");
1363  LogInfo("FastCalorimetry") << " RespCorrP (rsp) size " << rsp.size() << std::endl;
1364 
1365  if( rsp.size()%3 !=0 ) {
1366  LogInfo("FastCalorimetry")
1367  << " RespCorrP size is wrong -> no corrections applied !!!"
1368  << std::endl;
1369 
1370  p_knots.push_back(14000.);
1371  k_e.push_back (1.);
1372  k_h.push_back (1.);
1373  }
1374  else {
1375  for(unsigned i = 0; i < rsp.size(); i += 3) {
1376  LogInfo("FastCalorimetry") << "i = " << i/3 << " p = " << rsp [i]
1377  << " k_e(p) = " << rsp[i+1]
1378  << " k_e(p) = " << rsp[i+2] << std::endl;
1379 
1380  p_knots.push_back(rsp[i]);
1381  k_e.push_back (rsp[i+1]);
1382  k_h.push_back (rsp[i+2]);
1383  }
1384  }
1385 
1386 
1387  //FR
1388  edm::ParameterSet HCALparameters = fastCalo.getParameter<edm::ParameterSet>("HCAL");
1389  optionHDSim_ = HCALparameters.getParameter<int>("SimOption");
1390  hdGridSize_ = HCALparameters.getParameter<int>("GridSize");
1391  hdSimMethod_ = HCALparameters.getParameter<int>("SimMethod");
1392  //RF
1393 
1394  unfoldedMode_ = fastCalo.getUntrackedParameter<bool>("UnfoldedMode",false);
1395 }
1396 
1397 
1398 void CalorimetryManager::updateMap(uint32_t cellid,float energy,int id,std::map<uint32_t,std::vector<std::pair<int,float> > > & mymap)
1399 {
1400  // std::cout << " updateMap " << std::endl;
1401  std::map<unsigned,std::vector<std::pair<int,float> > >::iterator cellitr;
1402  cellitr = mymap.find(cellid);
1403  if(!unfoldedMode_) id=0;
1404  if( cellitr==mymap.end())
1405  {
1406  std::vector<std::pair<int,float> > myElement;
1407  myElement.push_back(std::pair<int,float> (id,energy));
1408  mymap[cellid]=myElement;
1409  }
1410  else
1411  {
1412  if(!unfoldedMode_)
1413  {
1414  cellitr->second[0].second+=energy;
1415  }
1416  else
1417  cellitr->second.push_back(std::pair<int,float>(id,energy));
1418  }
1419 }
1420 
1421 void CalorimetryManager::updateMap(int hi,float energy,int tid,std::vector<std::vector<std::pair<int,float> > > & mymap, std::vector<int>& firedCells)
1422 {
1423  // Standard case first : one entry per cell
1424  if(!unfoldedMode_)
1425  {
1426  // if new entry, update the list
1427  if(mymap[hi].size()==0)
1428  {
1429  firedCells.push_back(hi);
1430  mymap[hi].push_back(std::pair<int,float>(0,energy));
1431  }
1432  else
1433  mymap[hi][0].second+=energy;
1434  }
1435  else
1436  {
1437  // std::cout << "update map " << mymap[hi].size() << " " << hi << std::setw(8) << std::setprecision(6) << energy ;
1438  // std::cout << " " << mymap[hi][0].second << std::endl;
1439  // the minimal size is always 1 ; in this case, no push_back
1440  if(mymap[hi].size()==0)
1441  {
1442  // if(tid==0) std::cout << " Shit ! " << std::endl;
1443  firedCells.push_back(hi);
1444  }
1445 
1446  mymap[hi].push_back(std::pair<int,float>(tid,energy));
1447  }
1448 
1449 }
1450 
1451 
1453 
1454  int sizeP = p_knots.size();
1455 
1456  if(sizeP <= 1) {
1457  ecorr = 1.;
1458  hcorr = 1.;
1459  }
1460  else {
1461  int ip = -1;
1462  for (int i = 0; i < sizeP; i++) {
1463  if (p < p_knots[i]) { ip = i; break;}
1464  }
1465  if (ip == 0) {
1466  ecorr = k_e[0];
1467  hcorr = k_h[0];
1468  }
1469  else {
1470  if(ip == -1) {
1471  ecorr = k_e[sizeP-1];
1472  hcorr = k_h[sizeP-1];
1473  }
1474  else {
1475  double x1 = p_knots[ip-1];
1476  double x2 = p_knots[ip];
1477  double y1 = k_e[ip-1];
1478  double y2 = k_e[ip];
1479 
1480  if(x1 == x2) {
1481  // std::cout << " equal p_knots values!!! " << std::endl;
1482  }
1483 
1484  ecorr = (y1 + (y2 - y1) * (p - x1)/(x2 - x1));
1485 
1486  y1 = k_h[ip-1];
1487  y2 = k_h[ip];
1488  hcorr = (y1 + (y2 - y1) * (p - x1)/(x2 - x1));
1489 
1490  }
1491  }
1492  }
1493 
1494  if(debug_)
1495  LogDebug("FastCalorimetry") << " p, ecorr, hcorr = " << p << " "
1496  << ecorr << " " << hcorr << std::endl;
1497 
1498 }
1499 
1500 
1502 {
1503  unsigned size=firedCellsEB_.size();
1504  // float sum=0.;
1505  for(unsigned ic=0;ic<size;++ic)
1506  {
1507  int hi=firedCellsEB_[ic];
1508  if(!unfoldedMode_)
1509  {
1510  c.push_back(PCaloHit(EBDetId::unhashIndex(hi),EBMapping_[hi][0].second,0.,0));
1511  // std::cout << "Adding " << hi << " " << EBDetId::unhashIndex(hi) << " " ;
1512  // std::cout << EBMapping_[hi][0].second << " " << EBMapping_[hi][0].first << std::endl;
1513  }
1514  else
1515  {
1516  unsigned npart=EBMapping_[hi].size();
1517  for(unsigned ip=0;ip<npart;++ip)
1518  {
1519  c.push_back(PCaloHit(EBDetId::unhashIndex(hi),EBMapping_[hi][ip].second,0.,
1520  EBMapping_[hi][ip].first));
1521 
1522  }
1523  }
1524 
1525  // sum+=cellit->second;
1526  }
1527 
1528 // for(unsigned ic=0;ic<61200;++ic)
1529 // {
1530 // EBDetId myCell(EBDetId::unhashIndex(ic));
1531 // if(!myCell.null())
1532 // {
1533 // float total=0.;
1534 // for(unsigned id=0;id<EBMapping_[ic].size();++id)
1535 // total+=EBMapping_[ic][id].second;
1536 // if(EBMapping_[ic].size()>0)
1537 // std::cout << "Adding " << ic << " " << myCell << " " << std::setprecision(8) <<total << std::endl;
1538 // }
1539 // }
1540 
1541 
1542  // std::cout << " SUM : " << sum << std::endl;
1543  // std::cout << " Added " <<c.size() << " hits " <<std::endl;
1544 }
1545 
1546 
1548 {
1549  unsigned size=firedCellsEE_.size();
1550  // float sum=0.;
1551  for(unsigned ic=0;ic<size;++ic)
1552  {
1553  int hi=firedCellsEE_[ic];
1554  if(!unfoldedMode_)
1555  c.push_back(PCaloHit(EEDetId::unhashIndex(hi),EEMapping_[hi][0].second,0.,0));
1556  else
1557  {
1558  unsigned npart=EEMapping_[hi].size();
1559  for(unsigned ip=0;ip<npart;++ip)
1560  c.push_back(PCaloHit(EEDetId::unhashIndex(hi),EEMapping_[hi][ip].second,0.,
1561  EEMapping_[hi][ip].first));
1562  }
1563 
1564  // sum+=cellit->second;
1565  }
1566  // std::cout << " SUM : " << sum << std::endl;
1567  // std::cout << " Added " <<c.size() << " hits " <<std::endl;
1568 }
1569 
1571 {
1572  unsigned size=firedCellsHCAL_.size();
1573  // float sum=0.;
1574  for(unsigned ic=0;ic<size;++ic)
1575  {
1576  int hi=firedCellsHCAL_[ic];
1577  if(!unfoldedMode_)
1578  c.push_back(PCaloHit(theDetIds_[hi],HMapping_[hi][0].second,0.,0));
1579  else
1580  {
1581  unsigned npart=HMapping_[hi].size();
1582  for(unsigned ip=0;ip<npart;++ip)
1583  c.push_back(PCaloHit(theDetIds_[hi],HMapping_[hi][ip].second,0.,
1584  HMapping_[hi][ip].first));
1585  }
1586 
1587  // sum+=cellit->second;
1588  }
1589  // std::cout << " SUM : " << sum << std::endl;
1590  // std::cout << " Added " <<c.size() << " hits " <<std::endl;
1591 }
1592 
1594 {
1595  std::map<uint32_t,std::vector<std::pair< int,float> > >::const_iterator cellit;
1596  std::map<uint32_t,std::vector<std::pair <int,float> > >::const_iterator preshEnd=ESMapping_.end();
1597 
1598  for(cellit=ESMapping_.begin();cellit!=preshEnd;++cellit)
1599  {
1600  if(!unfoldedMode_)
1601  c.push_back(PCaloHit(cellit->first,cellit->second[0].second,0.,0));
1602  else
1603  {
1604  unsigned npart=cellit->second.size();
1605  for(unsigned ip=0;ip<npart;++ip)
1606  {
1607  c.push_back(PCaloHit(cellit->first,cellit->second[ip].second,0.,cellit->second[ip].first));
1608  }
1609  }
1610  }
1611 }
1612 
1613 
1614 // The main danger in this method is to screw up to relationships between particles
1615 // So, the muon FSimTracks created by FSimEvent.cc are simply to be updated
1617 {
1618  unsigned size=muons.size();
1619  for(unsigned i=0; i<size;++i)
1620  {
1621  int id=muons[i].trackId();
1622  if(abs(muons[i].type())!=13) continue;
1623  // identify the corresponding muon in the local collection
1624 
1625  std::vector<FSimTrack>::const_iterator itcheck=find_if(muonSimTracks.begin(),muonSimTracks.end(),FSimTrackEqual(id));
1626  if(itcheck!=muonSimTracks.end())
1627  {
1628  muons[i].setTkPosition(itcheck->trackerSurfacePosition());
1629  muons[i].setTkMomentum(itcheck->trackerSurfaceMomentum());
1630 // std::cout << " Found the SimTrack " << std::endl;
1631 // std::cout << *itcheck << std::endl;
1632 // std::cout << "SimTrack Id "<< id << " " << muons[i] << " " << std::endl;
1633  }
1634 // else
1635 // {
1636 // std::cout << " Calorimetery Manager : this should really not happen " << std::endl;
1637 // std::cout << " Was looking for " << id << " " << muons[i] << std::endl;
1638 // for(unsigned i=0;i<muonSimTracks.size();++i)
1639 // std::cout << muonSimTracks[i] << std::endl;
1640 // }
1641  }
1642 
1643 }
1644 
#define LogDebug(id)
void setSpotEnergy(double e)
Set the spot energy.
Definition: HcalHitMaker.h:28
double getHCALEnergyResponse(double e, int hit)
type
Definition: HCALResponse.h:22
T getParameter(std::string const &) const
EventNumber_t event() const
Definition: EventID.h:44
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
std::vector< double > rsp
std::vector< double > k_h
bool noEndVertex() const
no end vertex
RawParticle myElec
A few pointers to save time.
int hashedIndex() const
get a compact index for arrays
Definition: EBDetId.h:73
float charge() const
charge
Definition: FSimTrack.h:47
std::vector< PCaloHit > PCaloHitContainer
const ECALProperties * ecalProperties(int onEcal) const
ECAL properties.
Definition: Calorimeter.cc:73
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:129
int onLayer2() const
Definition: FSimTrack.h:96
GflashPiKShowerProfile * thePiKProfile
void bookByNumber(const std::string &name, int n1, int n2, int nx, float xmin, float xmax, int ny=0, float ymin=0., float ymax=0.)
Definition: Histos.cc:200
void updateMap(uint32_t cellid, float energy, int id, std::map< uint32_t, std::vector< std::pair< int, float > > > &mymap)
std::map< uint32_t, std::vector< std::pair< int, float > > > ESMapping_
double radLenIncm() const
Radiation length in cm.
GflashTrajectory * getHelix()
Definition: GflashShowino.h:30
MaterialEffects * theMuonEcalEffects
MaterialEffects * theMuonHcalEffects
const XYZTLorentzVector & momentum() const
Temporary (until move of SimTrack to Mathcore) - No! Actually very useful.
Definition: FSimTrack.h:168
void setCrackPadSurvivalProbability(double val)
Definition: EcalHitMaker.h:134
std::vector< std::vector< std::pair< int, float > > > EBMapping_
response responseHCAL(int mip, double energy, double eta, int partype)
GflashHadronShowerProfile * theProfile
const RawParticle & layer1Entrance() const
The particle at Preshower Layer 1.
Definition: FSimTrack.h:117
bool compute()
Compute the shower longitudinal and lateral development.
Definition: HDShower.cc:471
std::pair< double, double > response
Definition: HCALResponse.h:20
#define abs(x)
Definition: mlp_lapack.h:159
void reconstructECAL(const FSimTrack &track)
void setPreshower(PreshowerHitMaker *const myPresh)
set the preshower address
Definition: EMShower.cc:632
void loadFromPreshower(edm::PCaloHitContainer &c) const
#define NULL
Definition: scimark2.h:8
double npart
Definition: HydjetWrapper.h:45
const CaloSubdetectorGeometry * getHcalGeometry() const
Definition: Calorimeter.h:57
void readParameters(const edm::ParameterSet &fastCalo)
TRandom random
Definition: MVATrainer.cc:138
const PreshowerLayer1Properties * layer1Properties(int onLayer1) const
Preshower Layer1 properties.
Definition: Calorimeter.cc:100
void HDShowerSimulation(const FSimTrack &myTrack)
Hadronic Shower Simulation.
static EEDetId unhashIndex(int hi)
Definition: EEDetId.cc:115
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
void MuonMipSimulation(const FSimTrack &myTrack)
std::vector< std::vector< std::pair< int, float > > > EEMapping_
std::vector< double > p_knots
double charge(const std::vector< uint8_t > &Ampls)
void compute()
Compute the shower longitudinal and lateral development.
Definition: EMShower.cc:244
std::vector< int > firedCellsEB_
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)
std::vector< int > firedCellsHCAL_
int onEcal() const
Definition: FSimTrack.h:101
std::vector< FSimTrack > muonSimTracks
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
const RandomEngine * random
U second(std::pair< T, U > const &p)
const LandauFluctuationGenerator * aLandauGenerator
MuonBremsstrahlungSimulator * muonBremsstrahlungSimulator() const
Return the Muon Bremsstrahlung engine.
bool compute()
Compute the shower longitudinal and lateral development.
Definition: HFShower.cc:475
void setRadiusFactor(double r)
Definition: EcalHitMaker.h:130
GflashAntiProtonShowerProfile * theAntiProtonProfile
math::XYZVector XYZVector
int getmip()
Definition: HDShower.h:41
GflashProtonShowerProfile * theProtonProfile
const PreshowerLayer2Properties * layer2Properties(int onLayer2) const
Preshower Layer2 properties.
Definition: Calorimeter.cc:108
const T & max(const T &a, const T &b)
double getPathLengthAtShower()
Definition: GflashShowino.h:26
void book(const std::string &name, int nx, float xmin, float xmax, int ny=0, float ymin=0., float ymax=0.)
Book an histogram (1D or 2D)
Definition: Histos.cc:23
tuple result
Definition: query.py:137
void setTrackParameters(const XYZNormal &normal, double X0depthoffset, const FSimTrack &theTrack)
void loadFromEcalEndcap(edm::PCaloHitContainer &c) const
void setID(const int id)
Definition: RawParticle.cc:101
const XYZTLorentzVector & deltaMom() const
Returns the actual energy lost.
std::vector< DetId > theDetIds_
std::vector< double > radiusPreshowerCorrections_
std::vector< int > firedCellsEE_
const HCALProperties * hcalProperties(int onHcal) const
HCAL properties.
Definition: Calorimeter.cc:84
static Histos * instance()
Definition: Histos.cc:15
const RawParticle & ecalEntrance() const
The particle at ECAL entrance.
Definition: FSimTrack.h:123
int onVFcal() const
Definition: FSimTrack.h:111
const std::map< uint32_t, float > & getHits()
Definition: HcalHitMaker.h:37
std::vector< double > theTailIntervals_
bool first
Definition: L1TdeRCT.cc:79
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:39
std::vector< std::vector< std::pair< int, float > > > HMapping_
void loadMuonSimTracks(edm::SimTrackContainer &m) const
unsigned int nTracks() const
Number of tracks.
Definition: FSimEvent.cc:48
HCALResponse * myHDResponse_
CaloGeometryHelper * myCalorimeter_
const std::map< unsigned, float > & getHits()
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:39
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:284
Definition: DetId.h:20
void put(const std::string &file, std::string name="")
Write one or all histogram(s) in a file.
Definition: Histos.cc:69
void loadFromHcal(edm::PCaloHitContainer &c) const
#define M_PI
Definition: BFit3D.cc:3
void updateState(ParticlePropagator &myTrack, double radlen)
Compute the material effect (calls the sub class)
int hashedIndex() const
Definition: EEDetId.h:177
void EMShowerSimulation(const FSimTrack &myTrack)
Log< T >::type log(const T &t)
Definition: Log.h:22
GammaFunctionGenerator * aGammaGenerator
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:31
CLHEP::Hep3Vector Gflash3Vector
Definition: Gflash3Vector.h:6
std::vector< double > mipValues_
static std::vector< std::pair< int, float > > myZero_
static EBDetId unhashIndex(int hi)
get a DetId from a compact index for arrays
Definition: EBDetId.cc:12
bool null() const
is this a null id ?
Definition: DetId.h:47
CalorimeterNumber getCalorimeterNumber(const Gflash3Vector position)
void setTkPosition(const math::XYZVectorD &pos)
Definition: SimTrack.h:40
double flatShoot(double xmin=0.0, double xmax=1.0) const
Definition: RandomEngine.h:30
double getPathLengthOnEcal()
Definition: GflashShowino.h:25
math::XYZVector XYZPoint
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:40
tuple muons
Definition: patZpeak.py:38
std::vector< unsigned int > evtsToDebug_
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:86
bool preshowerPresent() const
void setPulledPadSurvivalProbability(double val)
Definition: EcalHitMaker.h:132
const RawParticle & layer2Entrance() const
The particle at Preshower Layer 2.
Definition: FSimTrack.h:120
int onLayer1() const
Definition: FSimTrack.h:91
std::vector< double > k_e
std::vector< GflashHit > & getGflashHitList()
int onHcal() const
Definition: FSimTrack.h:106
void setTkMomentum(const math::XYZTLorentzVectorD &mom)
Definition: SimTrack.h:42
const RawParticle & hcalEntrance() const
The particle at HCAL entrance.
Definition: FSimTrack.h:126
void setMipEnergy(double e1, double e2)
tuple cout
Definition: gather_cfg.py:41
void setPreshowerPresent(bool ps)
Definition: EcalHitMaker.h:137
double getHCALEnergyResolution(double e, int hit)
string s
Definition: asciidump.py:422
DetId getClosestCell(const XYZPoint &point, bool ecal, bool central) const
tuple status
Definition: ntuplemaker.py:245
EnergyLossSimulator * energyLossSimulator() const
Return the Energy Loss engine.
const std::map< uint32_t, float > & getHits()
not been done.
bool addHit(double r, double phi, unsigned layer=0)
add the hit in the HCAL in local coordinates
Definition: HcalHitMaker.cc:28
HSParameters * myHSParameters_
std::vector< SimTrack > SimTrackContainer
tuple size
Write out results.
void setVertex(const XYZTLorentzVector &vtx)
set the vertex
Definition: RawParticle.h:287
bool setDepth(double, bool inCm=false)
set the depth in X0 or Lambda0 units depending on showerType
int hashed_index() const
Definition: HcalDetId.cc:119
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
bool computeShower()
Definition: HDRShower.cc:46
void reconstructHCAL(const FSimTrack &myTrack)
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:642
void setGrid(EcalHitMaker *const myGrid)
set the grid address
Definition: EMShower.h:60
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:57
void loadFromEcalBarrel(edm::PCaloHitContainer &c) const