CMS 3D CMS Logo

FastHFShowerLibrary.cc
Go to the documentation of this file.
1 // File: FastHFShowerLibrary.cc
3 // Description: Shower library for Very forward hadron calorimeter
5 
21 
22 #include "Randomize.hh"
23 #include "CLHEP/Units/GlobalSystemOfUnits.h"
24 #include "CLHEP/Units/GlobalPhysicalConstants.h"
25 
26 // Geant4 headers
27 #include "G4ParticleDefinition.hh"
28 #include "G4DynamicParticle.hh"
29 #include "G4DecayPhysics.hh"
30 #include "G4ParticleTable.hh"
31 #include "G4ParticleTypes.hh"
32 
33 // STL headers
34 #include <vector>
35 #include <iostream>
36 #include <mutex>
37 
38 //#define DebugLog
39 
40 static std::once_flag initializeOnce;
41 
43  : fast(p) {
44  edm::ParameterSet m_HS = p.getParameter<edm::ParameterSet>("HFShowerLibrary");
45  applyFidCut = m_HS.getParameter<bool>("ApplyFiducialCut");
46 }
47 
49 
50  edm::LogInfo("FastCalorimetry") << "initHFShowerLibrary::initialization";
51 
53  iSetup.get<IdealGeometryRecord>().get(cpv);
54 
56  iSetup.get<HcalSimNumberingRecord>().get(hdc);
57  hcalConstants = (HcalDDDSimConstants*)(&(*hdc));
58 
59  std::string name = "HcalHits";
61  hfshower.reset(new HFShowerLibrary(name,*cpv,fast));
62 
63  //only one thread can be allowed to setup the G4 physics table.
64  std::call_once(initializeOnce,[]() {
65  // Geant4 particles
66  G4DecayPhysics decays;
67  decays.ConstructParticle();
68  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
69  partTable->SetReadiness();
70  });
71  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
72  hfshower->initRun(partTable, hcalConstants); // init particle code
73 }
74 
76 
77 #ifdef DebugLog
78  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: recoHFShowerLibrary ";
79 #endif
80 
81  if(!myTrack.onVFcal()) {
82 #ifdef DebugLog
83  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: we should not be here ";
84 #endif
85  }
86 
87  hitMap.clear();
88  double eGen = 1000.*myTrack.vfcalEntrance().e(); // energy in [MeV]
89  double delZv = (myTrack.vfcalEntrance().vertex().Z()>0.0) ? 50.0 : -50.0;
90  G4ThreeVector vertex( 10.*myTrack.vfcalEntrance().vertex().X(),
91  10.*myTrack.vfcalEntrance().vertex().Y(),
92  10.*myTrack.vfcalEntrance().vertex().Z()+delZv); // in [mm]
93 
94  G4ThreeVector direction(myTrack.vfcalEntrance().Vect().X(),
95  myTrack.vfcalEntrance().Vect().Y(),
96  myTrack.vfcalEntrance().Vect().Z());
97 
98  bool ok;
99  double weight = 1.0; // rad. damage
100  int parCode = myTrack.type();
101  double tSlice = 0.1*vertex.mag()/29.98;
102 
103  std::vector<HFShowerLibrary::Hit> hits =
104  hfshower->fillHits(vertex,direction,parCode,eGen,ok,weight,tSlice,false);
105 
106  for (unsigned int i=0; i<hits.size(); ++i) {
107  G4ThreeVector pos = hits[i].position;
108  int depth = hits[i].depth;
109  double time = hits[i].time;
110  if (!applyFidCut || (HFFibreFiducial::PMTNumber(pos)>0) ) {
111 // if (!applyFidCut || (applyFidCut && HFFibreFiducial::PMTNumber(pos)>0)) {
112  int det = 5;
113  int lay = 1;
114  uint32_t id = 0;
115  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, depth, lay);
116  id = numberingScheme.getUnitID(tmp);
117 
118  CaloHitID current_id(id,time,myTrack.id());
119  std::map<CaloHitID,float>::iterator cellitr;
120  cellitr = hitMap.find(current_id);
121  if(cellitr==hitMap.end()) {
122  hitMap.insert(std::pair<CaloHitID,float>(current_id,1.0));
123  } else {
124  cellitr->second += 1.0;
125  }
126  } // end of isItinFidVolume check
127  } // end loop over hits
128 
129 }
130 
132 
133  HcalDetId hid(id);
134  if (hid.subdet() == HcalForward) {
135  int eta = hid.ieta();
136  int phi = hid.iphi();
137  if (hcalConstants->maxHFDepth(eta,phi) > 2) {
138  if (hid.depth() <= 2) {
139  int dep = (G4UniformRand() > 0.5) ? (2+hid.depth()) : hid.depth();
140  id = HcalDetId(HcalForward,eta,phi,dep).rawId();
141  }
142  }
143  }
144 }
T getParameter(std::string const &) const
int maxHFDepth(const int ieta, const int iphi) const
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:139
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:49
const edm::ParameterSet fast
void recoHFShowerLibrary(const FSimTrack &myTrack)
virtual uint32_t getUnitID(const HcalNumberingFromDDD::HcalID &id)
Definition: weight.py:1
void modifyDepth(uint32_t &id)
static int PMTNumber(const G4ThreeVector &pe_effect)
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
int depth() const
get the tower depth
Definition: HcalDetId.cc:108
HcalNumberingScheme numberingScheme
double Y() const
y of vertex
Definition: RawParticle.h:275
int ieta() const
get the cell ieta
Definition: HcalDetId.h:56
double Z() const
z of vertex
Definition: RawParticle.h:276
static std::once_flag initializeOnce
int onVFcal() const
Definition: FSimTrack.h:111
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
int iphi() const
get the cell iphi
Definition: HcalDetId.cc:103
FastHFShowerLibrary(edm::ParameterSet const &p)
const T & get() const
Definition: EventSetup.h:56
double X() const
x of vertex
Definition: RawParticle.h:274
std::vector< std::vector< double > > tmp
Definition: MVATrainer.cc:100
int type() const
particle type (HEP PDT convension)
Definition: CoreSimTrack.h:25
std::unique_ptr< HFShowerLibrary > hfshower
int id() const
the index in FBaseSimEvent and other vectors
Definition: FSimTrack.h:86
void const initHFShowerLibrary(const edm::EventSetup &)
std::map< CaloHitID, float > hitMap
HcalDDDSimConstants * hcalConstants