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