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 
22 #include "G4LogicalVolumeStore.hh"
23 #include "G4LogicalVolume.hh"
24 #include "G4Step.hh"
25 #include "G4Track.hh"
26 #include "G4ParticleTable.hh"
27 #include "G4VProcess.hh"
28 #include "G4Trap.hh"
29 
30 #include <iostream>
31 #include <fstream>
32 #include <iomanip>
33 
34 //#define EDM_ML_DEBUG
35 
36 HGCSD::HGCSD(G4String name, const DDCompactView & cpv,
37  const SensitiveDetectorCatalog & clg,
38  edm::ParameterSet const & p, const SimTrackManager* manager) :
39  CaloSD(name, cpv, clg, p, manager,
40  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
41  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
42  numberingScheme(0), mouseBite_(0), slopeMin_(0), levelT_(99) {
43 
45  eminHit = m_HGC.getParameter<double>("EminHit")*CLHEP::MeV;
46  rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
47  waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
48  angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
49  double waferSize = m_HGC.getUntrackedParameter<double>("WaferSize")*CLHEP::mm;
50  double mouseBite = m_HGC.getUntrackedParameter<double>("MouseBite")*CLHEP::mm;
51  mouseBiteCut_ = waferSize*tan(30.0*CLHEP::deg) - mouseBite;
52 
53  //this is defined in the hgcsens.xml
54  G4String myName(this->nameOfSD());
56  nameX = "HGCal";
57  if (myName.find("HitsEE")!=std::string::npos) {
59  nameX = "HGCalEESensitive";
60  } else if (myName.find("HitsHEfront")!=std::string::npos) {
62  nameX = "HGCalHESiliconSensitive";
63  } else if (myName.find("HitsHEback")!=std::string::npos) {
65  nameX = "HGCalHEScintillatorSensitive";
66  }
67 
68 #ifdef EDM_ML_DEBUG
69  edm::LogInfo("HGCSim")<< "**************************************************"
70  << "\n"
71  << "* *"
72  << "\n"
73  << "* Constructing a HGCSD with name " << name << "\n"
74  << "* *"
75  << "\n"
76  << "**************************************************";
77 #endif
78  edm::LogInfo("HGCSim") << "HGCSD:: Threshold for storing hits: " << eminHit
79  << " for " << nameX << " subdet " << myFwdSubdet_;
80  edm::LogInfo("HGCSim") << "Reject MosueBite Flag: " << rejectMB_
81  << " Size of wafer " << waferSize << " Mouse Bite "
82  << mouseBite << ":" << mouseBiteCut_ << " along "
83  << angles_.size() << " axes";
84 }
85 
87  if (numberingScheme) delete numberingScheme;
88  if (mouseBite_) delete mouseBite_;
89 }
90 
91 bool HGCSD::ProcessHits(G4Step * aStep, G4TouchableHistory * ) {
92 
93  NaNTrap( aStep ) ;
94 
95  if (aStep == NULL) {
96  return true;
97  } else {
98  double r = aStep->GetPreStepPoint()->GetPosition().perp();
99  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
100 #ifdef EDM_ML_DEBUG
101  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
102  bool notaMuon = (parCode == mupPDG || parCode == mumPDG ) ? false : true;
103  G4LogicalVolume* lv =
104  aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
105  edm::LogInfo("HGCSim") << "HGCSD: Hit from standard path from "
106  << lv->GetName() << " for Track "
107  << aStep->GetTrack()->GetTrackID() << " ("
108  << aStep->GetTrack()->GetDefinition()->GetParticleName()
109  << ":" << notaMuon << ") R = " << r << " Z = " << z
110  << " slope = " << r/z << ":" << slopeMin_;
111 #endif
112  // Apply fiducial cuts
113  if (r/z >= slopeMin_) {
114  if (getStepInfo(aStep)) {
115  if (hitExists() == false && edepositEM+edepositHAD>0.)
117  }
118  }
119  return true;
120  }
121 }
122 
123 double HGCSD::getEnergyDeposit(G4Step* aStep) {
124  double wt1 = getResponseWt(aStep->GetTrack());
125  double wt2 = aStep->GetTrack()->GetWeight();
126  double destep = wt1*(aStep->GetTotalEnergyDeposit());
127  if (wt2 > 0) destep *= wt2;
128  return destep;
129 }
130 
131 uint32_t HGCSD::setDetUnitId(G4Step * aStep) {
132 
133  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
134  const G4VTouchable* touch = preStepPoint->GetTouchable();
135 
136  //determine the exact position in global coordinates in the mass geometry
137  G4ThreeVector hitPoint = preStepPoint->GetPosition();
138  float globalZ=touch->GetTranslation(0).z();
139  int iz( globalZ>0 ? 1 : -1);
140 
141  //convert to local coordinates (=local to the current volume):
142  G4ThreeVector localpos = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
143 
144  //get the det unit id with
146 
147  int layer(0), module(0), cell(0);
149  layer = touch->GetReplicaNumber(0);
150  module = touch->GetReplicaNumber(1);
151  } else {
152  if (touch->GetHistoryDepth() == levelT_) {
153  layer = touch->GetReplicaNumber(0);
154  module = -1;
155  cell = -1;
156 #ifdef EDM_ML_DEBUG
157  edm::LogInfo("HGCSim") << "Depths: " << touch->GetHistoryDepth()
158  << " name " << touch->GetVolume(0)->GetName()
159  << " layer:module:cell " << layer << ":"
160  << module << ":" << cell << std::endl;
161 #endif
162  } else {
163  layer = touch->GetReplicaNumber(2);
164  module = touch->GetReplicaNumber(1);
165  cell = touch->GetReplicaNumber(0);
166  }
167 #ifdef EDM_ML_DEBUG
168  edm::LogInfo("HGCSim") << "Depths: " << touch->GetHistoryDepth() <<" name "
169  << touch->GetVolume(0)->GetName()
170  << ":" << touch->GetReplicaNumber(0) << " "
171  << touch->GetVolume(1)->GetName()
172  << ":" << touch->GetReplicaNumber(1) << " "
173  << touch->GetVolume(2)->GetName()
174  << ":" << touch->GetReplicaNumber(2) << " "
175  << " layer:module:cell " << layer << ":" << module
176  << ":" << cell <<" Material " << mat->GetName()<<":"
177  << aStep->GetPreStepPoint()->GetMaterial()->GetRadlen()
178  << std::endl;
179 #endif
180  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.) return 0;
181  }
182 
183  uint32_t id = setDetUnitId (subdet, layer, module, cell, iz, localpos);
184  if (rejectMB_ && m_mode != HGCalGeometryMode::Square && id != 0) {
185  int det, z, lay, wafer, type, ic;
186  HGCalTestNumbering::unpackHexagonIndex(id, det, z, lay, wafer, type, ic);
187 #ifdef EDM_ML_DEBUG
188  edm::LogInfo("HGCSim") << "ID " << std::hex << id << std::dec << " Decode "
189  << det << ":" << z << ":" << lay << ":" << wafer
190  << ":" << type << ":" << ic << std::endl;
191 #endif
192  if (mouseBite_->exclude(hitPoint, z, wafer)) id = 0;
193  }
194  return id;
195 }
196 
197 void HGCSD::update(const BeginOfJob * job) {
198 
199  const edm::EventSetup* es = (*job)();
201  es->get<IdealGeometryRecord>().get(nameX,hdc);
202  if (hdc.isValid()) {
203  const HGCalDDDConstants* hgcons = hdc.product();
204  m_mode = hgcons->geomMode();
205  slopeMin_ = hgcons->minSlope();
206  levelT_ = hgcons->levelTop();
209  } else {
210  edm::LogError("HGCSim") << "HCalSD : Cannot find HGCalDDDConstants for "
211  << nameX;
212  throw cms::Exception("Unknown", "HGCSD") << "Cannot find HGCalDDDConstants for " << nameX << "\n";
213  }
214 #ifdef EDM_ML_DEBUG
215  edm::LogInfo("HGCSim") << "HGCSD::Initialized with mode " << m_mode
216  << " Slope cut " << slopeMin_ << " top Level "
217  << levelT_ << std::endl;
218 #endif
219 }
220 
222  G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable();
223  G4String particleName;
224  mumPDG = theParticleTable->FindParticle(particleName="mu-")->GetPDGEncoding();
225  mupPDG = theParticleTable->FindParticle(particleName="mu+")->GetPDGEncoding();
226 #ifdef EDM_ML_DEBUG
227  edm::LogInfo("HGCSim") << "HGCSD: Particle code for mu- = " << mumPDG
228  << " for mu+ = " << mupPDG;
229 #endif
230 }
231 
232 bool HGCSD::filterHit(CaloG4Hit* aHit, double time) {
233  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit));
234 }
235 
236 uint32_t HGCSD::setDetUnitId (ForwardSubdetector &subdet, int layer, int module,
237  int cell, int iz, G4ThreeVector &pos) {
238  uint32_t id = numberingScheme ?
239  numberingScheme->getUnitID(subdet, layer, module, cell, iz, pos) : 0;
240  return id;
241 }
242 
243 int HGCSD::setTrackID (G4Step* aStep) {
244  theTrack = aStep->GetTrack();
245 
246  double etrack = preStepPoint->GetKineticEnergy();
247  TrackInformation * trkInfo = (TrackInformation *)(theTrack->GetUserInformation());
248  int primaryID = trkInfo->getIDonCaloSurface();
249  if (primaryID == 0) {
250 #ifdef EDM_ML_DEBUG
251  edm::LogInfo("HGCSim") << "HGCSD: Problem with primaryID **** set by "
252  << "force to TkID **** " <<theTrack->GetTrackID();
253 #endif
254  primaryID = theTrack->GetTrackID();
255  }
256 
257  if (primaryID != previousID.trackID())
258  resetForNewPrimary(preStepPoint->GetPosition(), etrack);
259 
260  return primaryID;
261 }
float edepositEM
Definition: CaloSD.h:121
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
virtual ~HGCSD()
Definition: HGCSD.cc:86
HGCSD(G4String, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCSD.cc:36
HGCalGeometryMode::GeometryMode m_mode
Definition: HGCSD.h:52
std::string nameOfSD()
int levelT_
Definition: HGCSD.h:59
HGCNumberingScheme * numberingScheme
Definition: HGCSD.h:53
int getIDonCaloSurface() const
double mouseBiteCut_
Definition: HGCSD.h:61
Definition: CaloSD.h:42
virtual bool ProcessHits(G4Step *, G4TouchableHistory *)
Definition: HGCSD.cc:91
#define NULL
Definition: scimark2.h:8
std::vector< double > angles_
Definition: HGCSD.h:62
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:197
ForwardSubdetector
type of data representation of DDCompactView
Definition: DDCompactView.h:90
bool waferRot_
Definition: HGCSD.h:60
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
G4int mupPDG
Definition: HGCSD.h:55
double eminHit
Definition: HGCSD.h:56
const double MeV
double slopeMin_
Definition: HGCSD.h:58
void resetForNewPrimary(const G4ThreeVector &, double)
Definition: CaloSD.cc:454
virtual void initRun()
Definition: HGCSD.cc:221
float edepositHAD
Definition: CaloSD.h:121
int trackID() const
Definition: CaloHitID.h:25
double minSlope() const
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HGCalGeometryMode::GeometryMode geomMode() const
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
bool exclude(G4ThreeVector &point, int zside, int wafer)
Definition: HGCMouseBite.cc:25
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
virtual double getEnergyDeposit(G4Step *)
Definition: HGCSD.cc:123
virtual bool filterHit(CaloG4Hit *, double)
Definition: HGCSD.cc:232
bool rejectMB_
Definition: HGCSD.h:60
std::string nameX
Definition: HGCSD.h:50
const T & get() const
Definition: EventSetup.h:56
int setTrackID(G4Step *step)
Definition: HGCSD.cc:243
int levelTop() const
double getResponseWt(G4Track *)
Definition: CaloSD.cc:607
HLT enums.
G4bool hitExists()
Definition: CaloSD.cc:310
ForwardSubdetector myFwdSubdet_
Definition: HGCSD.h:57
Definition: vlib.h:208
virtual uint32_t setDetUnitId(G4Step *step)
Definition: HGCSD.cc:131
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
CaloG4Hit * createNewHit()
Definition: CaloSD.cc:363
G4int mumPDG
Definition: HGCSD.h:55
double getEnergyDeposit() const
Definition: CaloG4Hit.h:81
HGCMouseBite * mouseBite_
Definition: HGCSD.h:54