CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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) {
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  LogDebug("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 #ifdef DebugLog
85  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
86  bool notaMuon = (parCode == mupPDG || parCode == mumPDG ) ? false : true;
87  G4LogicalVolume* lv =
88  aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
89  edm::LogInfo("HGCSim") << "HGCSD: Hit from standard path from "
90  << lv->GetName() << " for Track "
91  << aStep->GetTrack()->GetTrackID() << " ("
92  << aStep->GetTrack()->GetDefinition()->GetParticleName()
93  << ":" << notaMuon << ")";
94 #endif
95  if (getStepInfo(aStep)) {
96  if (hitExists() == false && edepositEM+edepositHAD>0.) currentHit = createNewHit();
97  }
98  return true;
99  }
100 }
101 
102 double HGCSD::getEnergyDeposit(G4Step* aStep) {
103  double destep = aStep->GetTotalEnergyDeposit();
104  return destep;
105 }
106 
107 uint32_t HGCSD::setDetUnitId(G4Step * aStep) {
108 
109  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
110  const G4VTouchable* touch = preStepPoint->GetTouchable();
111 
112  //determine the exact position in global coordinates in the mass geometry
113  G4ThreeVector hitPoint = preStepPoint->GetPosition();
114  float globalZ=touch->GetTranslation(0).z();
115  int iz( globalZ>0 ? 1 : -1);
116 
117  //convert to local coordinates (=local to the current volume):
118  G4ThreeVector localpos = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
119 
120  //get the det unit id with
122 
123  int layer(0), module(0), cell(0);
125  layer = touch->GetReplicaNumber(0);
126  module = touch->GetReplicaNumber(1);
127  } else {
128  layer = touch->GetReplicaNumber(2);
129  module = touch->GetReplicaNumber(1);
130  cell = touch->GetReplicaNumber(0);
131  }
132 
133  return setDetUnitId (subdet, layer, module, cell, iz, localpos);
134 }
135 
136 void HGCSD::update(const BeginOfJob * job) {
137 
138  const edm::EventSetup* es = (*job)();
140  es->get<IdealGeometryRecord>().get(nameX,hdc);
141  if (hdc.isValid()) {
142  m_mode = hdc->geomMode();
144  } else {
145  edm::LogError("HGCSim") << "HCalSD : Cannot find HGCalDDDConstants for "
146  << nameX;
147  throw cms::Exception("Unknown", "HGCSD") << "Cannot find HGCalDDDConstants for " << nameX << "\n";
148  }
149 }
150 
152  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
153  G4String particleName;
154  mumPDG = theParticleTable->FindParticle(particleName="mu-")->GetPDGEncoding();
155  mupPDG = theParticleTable->FindParticle(particleName="mu+")->GetPDGEncoding();
156 #ifdef DebugLog
157  LogDebug("HGCSim") << "HGCSD: Particle code for mu- = " << mumPDG
158  << " for mu+ = " << mupPDG;
159 #endif
160 }
161 
162 bool HGCSD::filterHit(CaloG4Hit* aHit, double time) {
163  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit));
164 }
165 
166 
167 //
168 uint32_t HGCSD::setDetUnitId (ForwardSubdetector &subdet, int layer, int module,
169  int cell, int iz, G4ThreeVector &pos) {
170  return (numberingScheme ? numberingScheme->getUnitID(subdet, layer, module, cell, iz, pos) : 0);
171 }
172 
173 //
174 int HGCSD::setTrackID (G4Step* aStep) {
175  theTrack = aStep->GetTrack();
176 
177  double etrack = preStepPoint->GetKineticEnergy();
178  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
179  int primaryID = trkInfo->getIDonCaloSurface();
180  if (primaryID == 0) {
181 #ifdef DebugLog
182  edm::LogInfo("HGCSim") << "HGCSD: Problem with primaryID **** set by "
183  << "force to TkID **** " <<theTrack->GetTrackID();
184 #endif
185  primaryID = theTrack->GetTrackID();
186  }
187 
188  if (primaryID != previousID.trackID())
189  resetForNewPrimary(preStepPoint->GetPosition(), etrack);
190 
191  return primaryID;
192 }
#define LogDebug(id)
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()
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:136
ForwardSubdetector
type of data representation of DDCompactView
Definition: DDCompactView.h:77
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
void resetForNewPrimary(const G4ThreeVector &, double)
Definition: CaloSD.cc:454
virtual void initRun()
Definition: HGCSD.cc:151
float edepositHAD
Definition: CaloSD.h:121
int trackID() const
Definition: CaloHitID.h:25
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:102
virtual bool filterHit(CaloG4Hit *, double)
Definition: HGCSD.cc:162
std::string nameX
Definition: HGCSD.h:49
const T & get() const
Definition: EventSetup.h:56
int setTrackID(G4Step *step)
Definition: HGCSD.cc:174
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:107
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:363
G4int mumPDG
Definition: HGCSD.h:53
double getEnergyDeposit() const
Definition: CaloG4Hit.h:81