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 
16 
17 #include "Randomize.hh"
18 #include "CLHEP/Units/GlobalSystemOfUnits.h"
19 #include "CLHEP/Units/GlobalPhysicalConstants.h"
20 
21 // Geant4 headers
22 #include "G4ParticleDefinition.hh"
23 #include "G4DynamicParticle.hh"
24 #include "G4DecayPhysics.hh"
25 #include "G4ParticleTable.hh"
26 #include "G4ParticleTypes.hh"
27 
28 // STL headers
29 #include <vector>
30 #include <iostream>
31 
32 //#define DebugLog
33 
35  : fast(p)
36 {
37  edm::ParameterSet m_HS = p.getParameter<edm::ParameterSet>("HFShowerLibrary");
38  applyFidCut = m_HS.getParameter<bool>("ApplyFiducialCut");
39 }
40 
42 
43  edm::LogInfo("FastCalorimetry") << "initHFShowerLibrary::initialization";
44 
46  iSetup.get<IdealGeometryRecord>().get(cpv);
47 
48  std::string name = "HcalHits";
49  hfshower.reset(new HFShowerLibrary(name,*cpv,fast));
50  numberingFromDDD.reset(new HcalNumberingFromDDD(name, *cpv));
51 
52  // Geant4 particles
53  G4DecayPhysics decays;
54  decays.ConstructParticle();
55  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
56  partTable->SetReadiness();
57 
58  hfshower->initRun(partTable); // init particle code
59 }
60 
62 
63 #ifdef DebugLog
64  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: recoHFShowerLibrary ";
65 #endif
66 
67  if(!myTrack.onVFcal()) {
68 #ifdef DebugLog
69  edm::LogInfo("FastCalorimetry") << "FastHFShowerLibrary: we should not be here ";
70 #endif
71  }
72 
73  hitMap.clear();
74  double eGen = 1000.*myTrack.vfcalEntrance().e(); // energy in [MeV]
75  double delZv = (myTrack.vfcalEntrance().vertex().Z()>0.0) ? 50.0 : -50.0;
76  G4ThreeVector vertex( 10.*myTrack.vfcalEntrance().vertex().X(),
77  10.*myTrack.vfcalEntrance().vertex().Y(),
78  10.*myTrack.vfcalEntrance().vertex().Z()+delZv); // in [mm]
79 
80  G4ThreeVector direction(myTrack.vfcalEntrance().Vect().X(),
81  myTrack.vfcalEntrance().Vect().Y(),
82  myTrack.vfcalEntrance().Vect().Z());
83 
84  bool ok;
85  double weight = 1.0; // rad. damage
86  int parCode = myTrack.type();
87  double tSlice = 0.1*vertex.mag()/29.98;
88 
89  std::vector<HFShowerLibrary::Hit> hits =
90  hfshower->fillHits(vertex,direction,parCode,eGen,ok,weight,false,tSlice);
91 
92  for (unsigned int i=0; i<hits.size(); ++i) {
93  G4ThreeVector pos = hits[i].position;
94  int depth = hits[i].depth;
95  double time = hits[i].time;
96  if (!applyFidCut || (HFFibreFiducial::PMTNumber(pos)>0) ) {
97 // if (!applyFidCut || (applyFidCut && HFFibreFiducial::PMTNumber(pos)>0)) {
98  int det = 5;
99  int lay = 1;
100  uint32_t id = 0;
101  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, pos, depth, lay);
102  id = numberingScheme.getUnitID(tmp);
103 
104  CaloHitID current_id(id,time,myTrack.id());
105  std::map<CaloHitID,float>::iterator cellitr;
106  cellitr = hitMap.find(current_id);
107  if(cellitr==hitMap.end()) {
108  hitMap.insert(std::pair<CaloHitID,float>(current_id,1.0));
109  } else {
110  cellitr->second += 1.0;
111  }
112  } // end of isItinFidVolume check
113  } // end loop over hits
114 
115 }
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:55
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
int weight
Definition: histoStyle.py:50
void const initHFShowerLibrary(const edm::EventSetup &)
std::map< CaloHitID, float > hitMap