CMS 3D CMS Logo

HGCalSD.cc
Go to the documentation of this file.
1 // File: HGCalSD.cc
3 // Description: Sensitive Detector class for High Granularity Calorimeter
5 
13 #include "G4LogicalVolumeStore.hh"
14 #include "G4LogicalVolume.hh"
15 #include "G4Step.hh"
16 #include "G4Track.hh"
17 #include "G4ParticleTable.hh"
18 #include "G4VProcess.hh"
19 #include "G4Trap.hh"
20 
21 #include <fstream>
22 #include <iomanip>
23 #include <iostream>
24 #include <memory>
25 
26 //#define EDM_ML_DEBUG
27 
29  const HGCalDDDConstants* hgc,
30  const SensitiveDetectorCatalog& clg,
31  edm::ParameterSet const& p,
32  const SimTrackManager* manager)
33  : CaloSD(name,
34  clg,
35  p,
36  manager,
37  static_cast<float>(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
38  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
39  hgcons_(hgc),
40  slopeMin_(0),
41  levelT1_(99),
42  levelT2_(99),
43  tan30deg_(std::tan(30.0 * CLHEP::deg)) {
44  numberingScheme_.reset(nullptr);
45  mouseBite_.reset(nullptr);
46 
47  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCSD");
48  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
49  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
50  distanceFromEdge_ = m_HGC.getParameter<double>("DistanceFromEdge");
51  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
52  rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
53  waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
54  cornerMinMask_ = m_HGC.getParameter<int>("CornerMinMask");
55  angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
56 
57  if (storeAllG4Hits_) {
58  setUseMap(false);
60  }
61 
62  //this is defined in the hgcsens.xml
63  G4String myName = name;
65  nameX_ = "HGCal";
66  if (myName.find("HitsEE") != std::string::npos) {
68  nameX_ = "HGCalEESensitive";
69  } else if (myName.find("HitsHEfront") != std::string::npos) {
71  nameX_ = "HGCalHESiliconSensitive";
72  }
73 
74 #ifdef EDM_ML_DEBUG
75  edm::LogVerbatim("HGCSim") << "**************************************************"
76  << "\n"
77  << "* *"
78  << "\n"
79  << "* Constructing a HGCalSD with name " << name << "\n"
80  << "* *"
81  << "\n"
82  << "**************************************************";
83 #endif
84  edm::LogVerbatim("HGCSim") << "HGCalSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
85  << " detector " << mydet_;
86  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
87  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
88  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
89  edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_ << " cuts along " << angles_.size()
90  << " axes: " << angles_[0] << ", " << angles_[1];
91 }
92 
93 double HGCalSD::getEnergyDeposit(const G4Step* aStep) {
94  double r = aStep->GetPreStepPoint()->GetPosition().perp();
95  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
96 #ifdef EDM_ML_DEBUG
97  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
98  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
99  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
100  edm::LogVerbatim("HGCSim") << "HGCalSD: Hit from standard path from " << lv->GetName() << " for Track "
101  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
102  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
103 #endif
104  // Apply fiducial cut
105  if (r < z * slopeMin_) {
106 #ifdef EDM_ML_DEBUG
107  edm::LogVerbatim("HGCSim") << "HGCalSD: Fiducial Volume cut";
108 #endif
109  return 0.0;
110  }
111 
112  double wt1 = getResponseWt(aStep->GetTrack());
113  double wt2 = aStep->GetTrack()->GetWeight();
114  double destep = weight_ * wt1 * (aStep->GetTotalEnergyDeposit());
115  if (wt2 > 0)
116  destep *= wt2;
117 #ifdef EDM_ML_DEBUG
118  edm::LogVerbatim("HGCSim") << "HGCalSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << " Total weight "
119  << weight_ * wt1 * wt2 << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
120 #endif
121  return destep;
122 }
123 
124 uint32_t HGCalSD::setDetUnitId(const G4Step* aStep) {
125  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
126  const G4VTouchable* touch = preStepPoint->GetTouchable();
127 
128 #ifdef EDM_ML_DEBUG
129  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_;
130  printDetectorLevels(touch);
131 #endif
132  //determine the exact position in global coordinates in the mass geometry
133  G4ThreeVector hitPoint = preStepPoint->GetPosition();
134  float globalZ = touch->GetTranslation(0).z();
135  int iz(globalZ > 0 ? 1 : -1);
136 
137  int layer(0), module(-1), cell(-1);
139  if (touch->GetHistoryDepth() > levelT2_) {
140  layer = touch->GetReplicaNumber(4);
141  module = touch->GetReplicaNumber(3);
142  cell = touch->GetReplicaNumber(1);
143  } else {
144  layer = touch->GetReplicaNumber(3);
145  module = touch->GetReplicaNumber(2);
146  }
147 #ifdef EDM_ML_DEBUG
148  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
149  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
150  << module << ":" << cell;
151  printDetectorLevels(touch);
152 #endif
153  } else if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
154  layer = touch->GetReplicaNumber(0);
155 #ifdef EDM_ML_DEBUG
156  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
157  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
158  << module << ":" << cell;
159 #endif
160  } else {
161  layer = touch->GetReplicaNumber(3);
162  module = touch->GetReplicaNumber(2);
163  cell = touch->GetReplicaNumber(1);
164 #ifdef EDM_ML_DEBUG
165  edm::LogVerbatim("HGCSim") << "DepthsInside: " << touch->GetHistoryDepth() << " name "
166  << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":" << module
167  << ":" << cell;
168 #endif
169  }
170 #ifdef EDM_ML_DEBUG
171  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
172  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
173  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
174  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
175  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
176  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
177  << touch->GetReplicaNumber(4) << " "
178  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
179  << mat->GetName() << ":" << mat->GetRadlen();
180 #endif
181  // The following statement should be examined later before elimination
182  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
183  return 0;
184 
185  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
186  if (rejectMB_ && id != 0) {
187  auto uv = HGCSiliconDetId(id).waferUV();
188 #ifdef EDM_ML_DEBUG
189  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
190 #endif
191  if (mouseBite_->exclude(hitPoint, iz, uv.first, uv.second)) {
192  id = 0;
193 #ifdef EDM_ML_DEBUG
194  edm::LogVerbatim("HGCSim") << "Rejected by mousebite cutoff *****";
195 #endif
196  }
197  }
198 #ifdef EDM_ML_DEBUG
199  if (id != 0)
200  edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
201 #endif
202  return id;
203 }
204 
205 void HGCalSD::update(const BeginOfJob* job) {
206  if (hgcons_ != nullptr) {
209  levelT1_ = hgcons_->levelTop(0);
210  levelT2_ = hgcons_->levelTop(1);
211  double waferSize = hgcons_->waferSize(false);
212  double mouseBite = hgcons_->mouseBite(false);
213  mouseBiteCut_ = waferSize * tan30deg_ - mouseBite;
214 #ifdef EDM_ML_DEBUG
215  edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
216  << " top Level " << levelT1_ << ":" << levelT2_ << " wafer " << waferSize << ":"
217  << mouseBite;
218 #endif
219 
220  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_);
221  if (rejectMB_)
222  mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
223  } else {
224  throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
225  }
226 }
227 
229 
230 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
231  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
232 }
233 
234 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
235  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
236  if (cornerMinMask_ > 2) {
237  if (hgcons_->maskCell(DetId(id), cornerMinMask_)) {
238  id = 0;
239  ignoreRejection();
240  }
241  }
243  ignoreRejection();
244  return id;
245 }
246 
247 bool HGCalSD::isItinFidVolume(const G4ThreeVector& pos) {
248  if (fiducialCut_) {
249  return (hgcons_->distFromEdgeHex(pos.x(), pos.y(), pos.z()) > distanceFromEdge_);
250  } else {
251  return true;
252  }
253 }
bool maskCell(const DetId &id, int corners) const
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double distanceFromEdge_
Definition: HGCalSD.h:48
int levelT2_
Definition: HGCalSD.h:50
void setNumberCheckedHits(int val)
Definition: CaloSD.h:122
double slopeMin_
Definition: HGCalSD.h:48
Definition: CaloSD.h:40
HGCalSD(const std::string &, const HGCalDDDConstants *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCalSD.cc:28
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCalSD.h:43
std::string nameX_
Definition: HGCalSD.h:46
double mouseBiteCut_
Definition: HGCalSD.h:49
std::pair< int, int > waferUV() const
HGCalGeometryMode::GeometryMode geomMode() const
void setUseMap(bool val)
Definition: CaloSD.h:111
void ignoreRejection()
Definition: CaloSD.h:108
constexpr std::array< uint8_t, layerIndexSize > layer
T getUntrackedParameter(std::string const &, T const &) const
void initRun() override
Definition: HGCalSD.cc:228
double getEnergyDeposit(const G4Step *) override
Definition: HGCalSD.cc:93
bool rejectMB_
Definition: HGCalSD.h:52
DetId::Detector mydet_
Definition: HGCalSD.h:45
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCalSD.cc:124
bool storeAllG4Hits_
Definition: HGCalSD.h:51
std::vector< double > angles_
Definition: HGCalSD.h:54
const HGCalDDDConstants * hgcons_
Definition: HGCalSD.h:42
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const double tan30deg_
Definition: HGCalSD.h:53
double mouseBite(bool reco) const
bool waferRot_
Definition: HGCalSD.h:52
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCalSD.cc:247
double tmaxHit
Definition: CaloSD.h:144
Definition: DetId.h:17
int cornerMinMask_
Definition: HGCalSD.h:50
double minSlope() const
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCalSD.h:47
double weight_
Definition: HGCalSD.h:49
bool fiducialCut_
Definition: HGCalSD.h:52
void printDetectorLevels(const G4VTouchable *) const
Definition: CaloSD.cc:1072
HLT enums.
bool filterHit(CaloG4Hit *, double) override
Definition: HGCalSD.cc:230
int levelT1_
Definition: HGCalSD.h:50
std::unique_ptr< HGCMouseBite > mouseBite_
Definition: HGCalSD.h:44
double eminHit_
Definition: HGCalSD.h:48
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:849
int levelTop(int ind=0) const
double waferSize(bool reco) const
double distFromEdgeHex(double x, double y, double z) const
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCalSD.cc:205