CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
FastHFShowerLibrary.cc
Go to the documentation of this file.
1 // File: FastHFShowerLibrary.cc
3 // Description: Shower library for Very forward hadron calorimeter
5 
19 
20 #include "Randomize.hh"
21 #include "CLHEP/Units/GlobalSystemOfUnits.h"
22 #include "CLHEP/Units/GlobalPhysicalConstants.h"
23 
24 // Geant4 headers
25 #include "G4ParticleDefinition.hh"
26 #include "G4DynamicParticle.hh"
27 #include "G4DecayPhysics.hh"
28 #include "G4ParticleTable.hh"
29 #include "G4ParticleTypes.hh"
30 
31 // STL headers
32 #include <vector>
33 #include <iostream>
34 
35 //#define DebugLog
36 
38  : fast(p) {
39  edm::ParameterSet m_HS = p.getParameter<edm::ParameterSet>("HFShowerLibrary");
40  applyFidCut = m_HS.getParameter<bool>("ApplyFiducialCut");
41 }
42 
44 
45  edm::LogInfo("FastCalorimetry") << "initHFShowerLibrary::initialization";
46 
48  iSetup.get<IdealGeometryRecord>().get(cpv);
49 
51  iSetup.get<HcalSimNumberingRecord>().get(hdc);
52  HcalDDDSimConstants *hcalConstants = (HcalDDDSimConstants*)(&(*hdc));
53 
54  std::string name = "HcalHits";
55  numberingFromDDD.reset(new HcalNumberingFromDDD(hcalConstants));
56  hfshower.reset(new HFShowerLibrary(name,*cpv,fast));
57 
58  // Geant4 particles
59  G4DecayPhysics decays;
60  decays.ConstructParticle();
61  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
62  partTable->SetReadiness();
63 
64  hfshower->initRun(partTable, hcalConstants); // init particle code
65 }
66 
68 
69 #ifdef DebugLog
70  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: recoHFShowerLibrary ";
71 #endif
72 
73  if(!myTrack.onVFcal()) {
74 #ifdef DebugLog
75  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: we should not be here ";
76 #endif
77  }
78 
79  hitMap.clear();
80  double eGen = 1000.*myTrack.vfcalEntrance().e(); // energy in [MeV]
81  double delZv = (myTrack.vfcalEntrance().vertex().Z()>0.0) ? 50.0 : -50.0;
82  G4ThreeVector vertex( 10.*myTrack.vfcalEntrance().vertex().X(),
83  10.*myTrack.vfcalEntrance().vertex().Y(),
84  10.*myTrack.vfcalEntrance().vertex().Z()+delZv); // in [mm]
85 
86  G4ThreeVector direction(myTrack.vfcalEntrance().Vect().X(),
87  myTrack.vfcalEntrance().Vect().Y(),
88  myTrack.vfcalEntrance().Vect().Z());
89 
90  bool ok;
91  double weight = 1.0; // rad. damage
92  int parCode = myTrack.type();
93  double tSlice = 0.1*vertex.mag()/29.98;
94 
95  std::vector<HFShowerLibrary::Hit> hits =
96  hfshower->fillHits(vertex,direction,parCode,eGen,ok,weight,false,tSlice);
97 
98  for (unsigned int i=0; i<hits.size(); ++i) {
99  G4ThreeVector pos = hits[i].position;
100  int depth = hits[i].depth;
101  double time = hits[i].time;
102  if (!applyFidCut || (HFFibreFiducial::PMTNumber(pos)>0) ) {
103 // if (!applyFidCut || (applyFidCut && HFFibreFiducial::PMTNumber(pos)>0)) {
104  int det = 5;
105  int lay = 1;
106  uint32_t id = 0;
107  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, depth, lay);
108  id = numberingScheme.getUnitID(tmp);
109 
110  CaloHitID current_id(id,time,myTrack.id());
111  std::map<CaloHitID,float>::iterator cellitr;
112  cellitr = hitMap.find(current_id);
113  if(cellitr==hitMap.end()) {
114  hitMap.insert(std::pair<CaloHitID,float>(current_id,1.0));
115  } else {
116  cellitr->second += 1.0;
117  }
118  } // end of isItinFidVolume check
119  } // end loop over hits
120 
121 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
const RawParticle & vfcalEntrance() const
The particle at VFCAL entrance.
Definition: FSimTrack.h:139
const edm::ParameterSet fast
void recoHFShowerLibrary(const FSimTrack &myTrack)
virtual uint32_t getUnitID(const HcalNumberingFromDDD::HcalID &id)
static int PMTNumber(const G4ThreeVector &pe_effect)
HcalNumberingScheme numberingScheme
double Y() const
y of vertex
Definition: RawParticle.h:275
double Z() const
z of vertex
Definition: RawParticle.h:276
int onVFcal() const
Definition: FSimTrack.h:111
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
const XYZTLorentzVector & vertex() const
the vertex fourvector
Definition: RawParticle.h:285
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