CMS 3D CMS Logo

HGCSD.cc
Go to the documentation of this file.
1 // File: HGCSD.cc
3 // Description: Sensitive Detector class for Combined Forward Calorimeter
5 
7 
20 #include "G4LogicalVolumeStore.hh"
21 #include "G4LogicalVolume.hh"
22 #include "G4Step.hh"
23 #include "G4Track.hh"
24 #include "G4ParticleTable.hh"
25 #include "G4VProcess.hh"
26 #include "G4Trap.hh"
27 
28 #include <iostream>
29 #include <fstream>
30 #include <iomanip>
31 
32 //#define DebugLog
33 
34 HGCSD::HGCSD(G4String name, const DDCompactView & cpv,
35  const SensitiveDetectorCatalog & clg,
36  edm::ParameterSet const & p, const SimTrackManager* manager) :
37  CaloSD(name, cpv, clg, p, manager,
38  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
39  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
40  numberingScheme(0), slopeMin_(0), levelT_(99) {
41 
43  eminHit = m_HGC.getParameter<double>("EminHit")*CLHEP::MeV;
44 
45  //this is defined in the hgcsens.xml
46  G4String myName(this->nameOfSD());
48  nameX = "HGCal";
49  if (myName.find("HitsEE")!=std::string::npos) {
51  nameX = "HGCalEESensitive";
52  } else if (myName.find("HitsHEfront")!=std::string::npos) {
54  nameX = "HGCalHESiliconSensitive";
55  } else if (myName.find("HitsHEback")!=std::string::npos) {
57  nameX = "HGCalHEScintillatorSensitive";
58  }
59 
60 #ifdef DebugLog
61  edm::LogInfo("HGCSim")<< "**************************************************"
62  << "\n"
63  << "* *"
64  << "\n"
65  << "* Constructing a HGCSD with name " << name << "\n"
66  << "* *"
67  << "\n"
68  << "**************************************************";
69 #endif
70  edm::LogInfo("HGCSim") << "HGCSD:: Threshold for storing hits: " << eminHit;
71 }
72 
74  if (numberingScheme) delete numberingScheme;
75 }
76 
77 bool HGCSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
78 
79  NaNTrap( aStep ) ;
80 
81  if (aStep == NULL) {
82  return true;
83  } else {
84  double r = aStep->GetPreStepPoint()->GetPosition().perp();
85  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
86 #ifdef DebugLog
87  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
88  bool notaMuon = (parCode == mupPDG || parCode == mumPDG ) ? false : true;
89  G4LogicalVolume* lv =
90  aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
91  edm::LogInfo("HGCSim") << "HGCSD: Hit from standard path from "
92  << lv->GetName() << " for Track "
93  << aStep->GetTrack()->GetTrackID() << " ("
94  << aStep->GetTrack()->GetDefinition()->GetParticleName()
95  << ":" << notaMuon << ") R = " << r << " Z = " << z
96  << " slope = " << r/z << ":" << slopeMin_;
97 #endif
98  // Apply fiducial cuts
99  if (r/z >= slopeMin_) {
100  if (getStepInfo(aStep)) {
101  if (hitExists() == false && edepositEM+edepositHAD>0.)
103  }
104  }
105  return true;
106  }
107 }
108 
109 double HGCSD::getEnergyDeposit(G4Step* aStep) {
110  double destep = aStep->GetTotalEnergyDeposit();
111  return destep;
112 }
113 
114 uint32_t HGCSD::setDetUnitId(G4Step * aStep) {
115 
116  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
117  const G4VTouchable* touch = preStepPoint->GetTouchable();
118 
119  //determine the exact position in global coordinates in the mass geometry
120  G4ThreeVector hitPoint = preStepPoint->GetPosition();
121  float globalZ=touch->GetTranslation(0).z();
122  int iz( globalZ>0 ? 1 : -1);
123 
124  //convert to local coordinates (=local to the current volume):
125  G4ThreeVector localpos = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
126 
127  //get the det unit id with
129 
130  int layer(0), module(0), cell(0);
132  layer = touch->GetReplicaNumber(0);
133  module = touch->GetReplicaNumber(1);
134  } else {
135  if (touch->GetHistoryDepth() == levelT_) {
136  layer = touch->GetReplicaNumber(0);
137  module = -1;
138  cell = -1;
139 #ifdef DebugLog
140  std::cout << "Depths: " << touch->GetHistoryDepth() << " name "
141  << touch->GetVolume(0)->GetLogicalVolume()->GetName()
142  << " layer:module:cell " << layer << ":" << module << ":"
143  << cell << std::endl;
144 #endif
145  } else {
146  layer = touch->GetReplicaNumber(2);
147  module = touch->GetReplicaNumber(1);
148  cell = touch->GetReplicaNumber(0);
149  }
150  }
151 
152  return setDetUnitId (subdet, layer, module, cell, iz, localpos);
153 }
154 
155 void HGCSD::update(const BeginOfJob * job) {
156 
157  const edm::EventSetup* es = (*job)();
159  es->get<IdealGeometryRecord>().get(nameX,hdc);
160  if (hdc.isValid()) {
161  m_mode = hdc->geomMode();
162  slopeMin_ = hdc->minSlope();
163  levelT_ = hdc->levelTop();
165  } else {
166  edm::LogError("HGCSim") << "HCalSD : Cannot find HGCalDDDConstants for "
167  << nameX;
168  throw cms::Exception("Unknown", "HGCSD") << "Cannot find HGCalDDDConstants for " << nameX << "\n";
169  }
170 #ifdef DebugLog
171  edm::LogInfo("HGCSim") << "HGCSD::Initialized with mode " << m_mode
172  << " Slope cut " << slopeMin_ << " top Level "
173  << levelT_ << std::endl;
174 #endif
175 }
176 
178  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
179  G4String particleName;
180  mumPDG = theParticleTable->FindParticle(particleName="mu-")->GetPDGEncoding();
181  mupPDG = theParticleTable->FindParticle(particleName="mu+")->GetPDGEncoding();
182 #ifdef DebugLog
183  edm::LogInfo("HGCSim") << "HGCSD: Particle code for mu- = " << mumPDG
184  << " for mu+ = " << mupPDG;
185 #endif
186 }
187 
188 bool HGCSD::filterHit(CaloG4Hit* aHit, double time) {
189  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit));
190 }
191 
192 uint32_t HGCSD::setDetUnitId (ForwardSubdetector &subdet, int layer, int module,
193  int cell, int iz, G4ThreeVector &pos) {
194  return (numberingScheme ? numberingScheme->getUnitID(subdet, layer, module, cell, iz, pos) : 0);
195 }
196 
197 int HGCSD::setTrackID (G4Step* aStep) {
198  theTrack = aStep->GetTrack();
199 
200  double etrack = preStepPoint->GetKineticEnergy();
201  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
202  int primaryID = trkInfo->getIDonCaloSurface();
203  if (primaryID == 0) {
204 #ifdef DebugLog
205  edm::LogInfo("HGCSim") << "HGCSD: Problem with primaryID **** set by "
206  << "force to TkID **** " <<theTrack->GetTrackID();
207 #endif
208  primaryID = theTrack->GetTrackID();
209  }
210 
211  if (primaryID != previousID.trackID())
212  resetForNewPrimary(preStepPoint->GetPosition(), etrack);
213 
214  return primaryID;
215 }
float edepositEM
Definition: CaloSD.h:121
T getParameter(std::string const &) const
virtual ~HGCSD()
Definition: HGCSD.cc:73
HGCSD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCSD.cc:34
std::string nameOfSD()
int levelT_
Definition: HGCSD.h:57
HGCNumberingScheme * numberingScheme
Definition: HGCSD.h:52
int getIDonCaloSurface() const
Definition: CaloSD.h:42
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: HGCSD.cc:77
#define NULL
Definition: scimark2.h:8
uint32_t getUnitID(ForwardSubdetector subdet, int layer, int module, int cell, int iz, const G4ThreeVector &pos)
assigns the det id to a hit
virtual void update(const BeginOfJob *)
This routine will be called when the appropriate signal arrives.
Definition: HGCSD.cc:155
ForwardSubdetector
type of data representation of DDCompactView
Definition: DDCompactView.h:90
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
G4int mupPDG
Definition: HGCSD.h:53
double eminHit
Definition: HGCSD.h:54
HGCalGeometryMode m_mode
Definition: HGCSD.h:51
const double MeV
double slopeMin_
Definition: HGCSD.h:56
void resetForNewPrimary(const G4ThreeVector &, double)
Definition: CaloSD.cc:454
virtual void initRun()
Definition: HGCSD.cc:177
float edepositHAD
Definition: CaloSD.h:121
int trackID() const
Definition: CaloHitID.h:25
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
CaloHitID previousID
Definition: CaloSD.h:117
CaloG4Hit * currentHit
Definition: CaloSD.h:128
void NaNTrap(G4Step *step)
G4Track * theTrack
Definition: CaloSD.h:118
virtual G4bool getStepInfo(G4Step *aStep)
Definition: CaloSD.cc:250
double tmaxHit
Definition: CaloSD.h:123
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
virtual double getEnergyDeposit(G4Step *)
Definition: HGCSD.cc:109
virtual bool filterHit(CaloG4Hit *, double)
Definition: HGCSD.cc:188
std::string nameX
Definition: HGCSD.h:49
const T & get() const
Definition: EventSetup.h:56
int setTrackID(G4Step *step)
Definition: HGCSD.cc:197
HLT enums.
G4bool hitExists()
Definition: CaloSD.cc:310
ForwardSubdetector myFwdSubdet_
Definition: HGCSD.h:55
Definition: vlib.h:208
virtual uint32_t setDetUnitId(G4Step *step)
Definition: HGCSD.cc:114
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:363
G4int mumPDG
Definition: HGCSD.h:53
double getEnergyDeposit() const
Definition: CaloG4Hit.h:81