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 
14 #include "G4LogicalVolumeStore.hh"
15 #include "G4LogicalVolume.hh"
16 #include "G4Step.hh"
17 #include "G4Track.hh"
18 #include "G4ParticleTable.hh"
19 #include "G4VProcess.hh"
20 #include "G4Trap.hh"
21 
22 #include <fstream>
23 #include <iomanip>
24 #include <iostream>
25 #include <memory>
26 
27 //#define EDM_ML_DEBUG
28 
29 using namespace angle_units::operators;
30 
32  const HGCalDDDConstants* hgc,
33  const SensitiveDetectorCatalog& clg,
34  edm::ParameterSet const& p,
35  const SimTrackManager* manager)
36  : CaloSD(name,
37  clg,
38  p,
39  manager,
40  static_cast<float>(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
41  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
42  hgcons_(hgc),
43  slopeMin_(0),
44  levelT1_(99),
45  levelT2_(99),
46  useSimWt_(0),
47  tan30deg_(std::tan(30.0 * CLHEP::deg)),
48  cos30deg_(std::cos(30.0 * CLHEP::deg)) {
49  numberingScheme_.reset(nullptr);
50  guardRing_.reset(nullptr);
51  guardRingPartial_.reset(nullptr);
52  mouseBite_.reset(nullptr);
53 
54  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCSD");
55  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
56  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
57  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
58  rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
59  waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
60  cornerMinMask_ = m_HGC.getParameter<int>("CornerMinMask");
61  angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
62  missingFile_ = m_HGC.getUntrackedParameter<std::string>("MissingWaferFile");
63  checkID_ = m_HGC.getUntrackedParameter<bool>("CheckID");
64  verbose_ = m_HGC.getUntrackedParameter<int>("Verbosity");
65 
66  if (storeAllG4Hits_) {
67  setUseMap(false);
69  }
70 
71  //this is defined in the hgcsens.xml
72  G4String myName = name;
74  nameX_ = "HGCal";
75  if (myName.find("HitsEE") != std::string::npos) {
77  nameX_ = "HGCalEESensitive";
78  } else if (myName.find("HitsHEfront") != std::string::npos) {
80  nameX_ = "HGCalHESiliconSensitive";
81  }
82 
83 #ifdef EDM_ML_DEBUG
84  edm::LogVerbatim("HGCSim") << "**************************************************"
85  << "\n"
86  << "* *"
87  << "\n"
88  << "* Constructing a HGCalSD with name " << name << "\n"
89  << "* *"
90  << "\n"
91  << "**************************************************";
92 #endif
93  edm::LogVerbatim("HGCSim") << "HGCalSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
94  << " detector " << mydet_;
95  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
96  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
97  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
98  edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_ << " cuts along " << angles_.size()
99  << " axes: " << angles_[0] << ", " << angles_[1];
100 }
101 
102 double HGCalSD::getEnergyDeposit(const G4Step* aStep) {
103  double r = aStep->GetPreStepPoint()->GetPosition().perp();
104  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
105 #ifdef EDM_ML_DEBUG
106  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
107  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
108  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
109  edm::LogVerbatim("HGCSim") << "HGCalSD: Hit from standard path from " << lv->GetName() << " for Track "
110  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
111  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
112 #endif
113  // Apply fiducial cut
114  if (r < z * slopeMin_) {
115 #ifdef EDM_ML_DEBUG
116  edm::LogVerbatim("HGCSim") << "HGCalSD: Fiducial Volume cut";
117 #endif
118  return 0.0;
119  }
120 
121  double wt1 = getResponseWt(aStep->GetTrack());
122  double wt2 = aStep->GetTrack()->GetWeight();
123  double destep = weight_ * wt1 * (aStep->GetTotalEnergyDeposit());
124  if (wt2 > 0)
125  destep *= wt2;
126 #ifdef EDM_ML_DEBUG
127  edm::LogVerbatim("HGCSim") << "HGCalSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << " Total weight "
128  << weight_ * wt1 * wt2 << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
129 #endif
130  return destep;
131 }
132 
133 uint32_t HGCalSD::setDetUnitId(const G4Step* aStep) {
134  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
135  const G4VTouchable* touch = preStepPoint->GetTouchable();
136 
137 #ifdef EDM_ML_DEBUG
138  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_;
139  printDetectorLevels(touch);
140 #endif
141  //determine the exact position in global coordinates in the mass geometry
142  G4ThreeVector hitPoint = preStepPoint->GetPosition();
143  float globalZ = touch->GetTranslation(0).z();
144  int iz(globalZ > 0 ? 1 : -1);
145 
146  int layer(0), moduleLev(-1), cell(-1);
147  if (useSimWt_ > 0) {
148  layer = touch->GetReplicaNumber(2);
149  moduleLev = 1;
150  } else if (touch->GetHistoryDepth() > levelT2_) {
151  layer = touch->GetReplicaNumber(4);
152  cell = touch->GetReplicaNumber(1);
153  moduleLev = 3;
154  } else {
155  layer = touch->GetReplicaNumber(3);
156  moduleLev = 2;
157  }
158  int module = touch->GetReplicaNumber(moduleLev);
159  if (verbose_ && (cell == -1))
160  edm::LogVerbatim("HGCSim") << "Top " << touch->GetVolume(0)->GetName() << " Module "
161  << touch->GetVolume(moduleLev)->GetName();
162 #ifdef EDM_ML_DEBUG
163  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << ":"
164  << useSimWt_ << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell "
165  << layer << ":" << moduleLev << ":" << module << ":" << cell;
166  printDetectorLevels(touch);
167  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
168  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
169  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
170  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
171  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
172  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
173  << touch->GetReplicaNumber(4) << " "
174  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
175  << mat->GetName() << ":" << mat->GetRadlen();
176 #endif
177  // The following statement should be examined later before elimination
178  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
179  return 0;
180 
181  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
182  if ((rejectMB_ || fiducialCut_) && id != 0) {
183  auto uv = HGCSiliconDetId(id).waferUV();
184 #ifdef EDM_ML_DEBUG
185  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
186 #endif
187  G4ThreeVector local = (touch->GetHistory()->GetTransform(moduleLev).TransformPoint(hitPoint));
188  if (fiducialCut_) {
189  int layertype = hgcons_->layerType(layer);
190  int frontBack = HGCalTypes::layerFrontBack(layertype);
191  if (guardRing_->exclude(local, iz, frontBack, layer, uv.first, uv.second) ||
192  guardRingPartial_->exclude(local, iz, frontBack, layer, uv.first, uv.second)) {
193  id = 0;
194 #ifdef EDM_ML_DEBUG
195  edm::LogVerbatim("HGCSim") << "Rejected by GuardRing cutoff *****";
196 #endif
197  }
198  }
199  if ((rejectMB_) && (mouseBite_->exclude(local, iz, layer, uv.first, uv.second))) {
200  id = 0;
201 #ifdef EDM_ML_DEBUG
202  edm::LogVerbatim("HGCSim") << "Rejected by MouseBite cutoff *****";
203 #endif
204  }
205  }
206 #ifdef EDM_ML_DEBUG
207  if (id != 0)
208  edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
209 #endif
210  if ((id != 0) && checkID_) {
211  HGCSiliconDetId hid1(id);
212  bool cshift = (hgcons_->cassetteShiftSilicon(hid1.zside(), hid1.layer(), hid1.waferU(), hid1.waferV()));
213  std::string_view pid = (cshift ? "HGCSim" : "HGCalSim");
214  bool debug = (verbose_ > 0) ? true : false;
215  auto xy = hgcons_->locateCell(hid1, debug);
216  double xx = (hid1.zside() > 0) ? xy.first : -xy.first;
217  double dx = xx - (hitPoint.x() / CLHEP::cm);
218  double dy = xy.second - (hitPoint.y() / CLHEP::cm);
219  double diff = (dx * dx + dy * dy);
220  constexpr double tol = 2.0 * 2.0;
221  bool valid1 = hgcons_->isValidHex8(hid1.layer(), hid1.waferU(), hid1.waferV(), hid1.cellU(), hid1.cellV(), true);
222  if ((diff > tol) || (!valid1))
223  pid = "HGCalError";
224  auto partn = hgcons_->waferTypeRotation(hid1.layer(), hid1.waferU(), hid1.waferV(), false, false);
225  int indx = HGCalWaferIndex::waferIndex(layer, hid1.waferU(), hid1.waferV());
226  edm::LogVerbatim(pid) << "CheckID " << HGCSiliconDetId(id) << " Layer:Module:Cell:ModuleLev " << layer << ":"
227  << module << ":" << cell << ":" << moduleLev << " SimWt:history " << useSimWt_ << ":"
228  << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << " input position: ("
229  << hitPoint.x() / CLHEP::cm << ", " << hitPoint.y() / CLHEP::cm << ":"
230  << convertRadToDeg(std::atan2(hitPoint.y(), hitPoint.x())) << "); position from ID (" << xx
231  << ", " << xy.second << ") distance " << dx << ":" << dy << ":" << std::sqrt(diff)
232  << " Valid " << valid1 << " Wafer type|rotation " << partn.first << ":" << partn.second
233  << " Part:Orient:Cassette " << std::get<1>(hgcons_->waferFileInfo(indx)) << ":"
234  << std::get<2>(hgcons_->waferFileInfo(indx)) << ":"
235  << std::get<3>(hgcons_->waferFileInfo(indx)) << " CassetteShift " << cshift;
236  if ((diff > tol) || (!valid1)) {
237  printDetectorLevels(touch);
238  hgcons_->locateCell(hid1, true);
239  }
240  }
241  return id;
242 }
243 
244 void HGCalSD::update(const BeginOfJob* job) {
245  if (hgcons_ != nullptr) {
248  levelT1_ = hgcons_->levelTop(0);
249  levelT2_ = hgcons_->levelTop(1);
251  int useOffset = hgcons_->getParameter()->useOffset_;
252  double waferSize = hgcons_->waferSize(false);
253  double mouseBite = hgcons_->mouseBite(false);
254  double guardRingOffset = hgcons_->guardRingOffset(false);
255  double sensorSizeOffset = hgcons_->sensorSizeOffset(false);
256  if (useOffset > 0) {
257  rejectMB_ = true;
258  fiducialCut_ = true;
259  }
260  double mouseBiteNew = (fiducialCut_) ? (mouseBite + guardRingOffset + sensorSizeOffset / cos30deg_) : mouseBite;
261  mouseBiteCut_ = waferSize * tan30deg_ - mouseBiteNew;
262 #ifdef EDM_ML_DEBUG
263  edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
264  << " top Level " << levelT1_ << ":" << levelT2_ << " useSimWt " << useSimWt_ << " wafer "
265  << waferSize << ":" << mouseBite << ":" << guardRingOffset << ":" << sensorSizeOffset
266  << ":" << mouseBiteNew << ":" << mouseBiteCut_ << " useOffset " << useOffset;
267 #endif
268 
269  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_, missingFile_);
270  if (rejectMB_)
271  mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
272  if (fiducialCut_) {
273  guardRing_ = std::make_unique<HGCGuardRing>(*hgcons_);
274  guardRingPartial_ = std::make_unique<HGCGuardRingPartial>(*hgcons_);
275  }
276  } else {
277  throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
278  }
279 }
280 
282 
283 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
284  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
285 }
286 
287 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
288  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
289  if (cornerMinMask_ > 2) {
290  if (hgcons_->maskCell(DetId(id), cornerMinMask_)) {
291  id = 0;
292  ignoreRejection();
293  }
294  }
295  if (hgcons_->waferHexagon8File() || (id == 0))
296  ignoreRejection();
297  return id;
298 }
bool maskCell(const DetId &id, int corners) const
Log< level::Info, true > LogVerbatim
int verbose_
Definition: HGCalSD.h:56
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double distanceFromEdge_
Definition: HGCalSD.h:51
int levelT2_
Definition: HGCalSD.h:53
void setNumberCheckedHits(int val)
Definition: CaloSD.h:122
double slopeMin_
Definition: HGCalSD.h:51
const HGCalParameters * getParameter() const
Definition: CaloSD.h:40
HGCalSD(const std::string &, const HGCalDDDConstants *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCalSD.cc:31
constexpr NumType convertRadToDeg(NumType radians)
Definition: angle_units.h:21
int cellU() const
get the cell #&#39;s in u,v or in x,y
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCalSD.h:44
std::pair< int, int > waferTypeRotation(int layer, int waferU, int waferV, bool fromFile, bool debug) const
std::string nameX_
Definition: HGCalSD.h:49
double mouseBiteCut_
Definition: HGCalSD.h:52
std::pair< int, int > waferUV() const
HGCalGeometryMode::GeometryMode geomMode() const
bool isValidHex8(int lay, int waferU, int waferV, bool fullAndPart) const
void setUseMap(bool val)
Definition: CaloSD.h:111
void ignoreRejection()
Definition: CaloSD.h:108
T getUntrackedParameter(std::string const &, T const &) const
int layer() const
get the layer #
void initRun() override
Definition: HGCalSD.cc:281
double getEnergyDeposit(const G4Step *) override
Definition: HGCalSD.cc:102
bool rejectMB_
Definition: HGCalSD.h:55
DetId::Detector mydet_
Definition: HGCalSD.h:48
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCalSD.cc:133
std::tuple< int, int, int, int > waferFileInfo(unsigned int kk) const
T sqrt(T t)
Definition: SSEVec.h:19
const double cos30deg_
Definition: HGCalSD.h:57
std::string missingFile_
Definition: HGCalSD.h:59
bool storeAllG4Hits_
Definition: HGCalSD.h:54
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
std::vector< double > angles_
Definition: HGCalSD.h:58
const HGCalDDDConstants * hgcons_
Definition: HGCalSD.h:43
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int waferU() const
const double tan30deg_
Definition: HGCalSD.h:57
double mouseBite(bool reco) const
bool waferRot_
Definition: HGCalSD.h:55
int32_t waferIndex(int32_t layer, int32_t waferU, int32_t waferV, bool old=false)
double tmaxHit
Definition: CaloSD.h:144
bool checkID_
Definition: HGCalSD.h:55
Definition: DetId.h:17
int cornerMinMask_
Definition: HGCalSD.h:53
#define debug
Definition: HDRShower.cc:19
double minSlope() const
int waferV() const
int zside() const
get the z-side of the cell (1/-1)
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCalSD.h:50
std::pair< float, float > locateCell(int cell, int lay, int type, bool reco) const
bool waferHexagon8File() const
int cellV() const
double weight_
Definition: HGCalSD.h:52
bool fiducialCut_
Definition: HGCalSD.h:55
void printDetectorLevels(const G4VTouchable *) const
Definition: CaloSD.cc:1083
std::unique_ptr< HGCGuardRingPartial > guardRingPartial_
Definition: HGCalSD.h:46
HLT enums.
bool filterHit(CaloG4Hit *, double) override
Definition: HGCalSD.cc:283
bool cassetteShiftSilicon(int zside, int layer, int waferU, int waferV) const
double guardRingOffset(bool reco) const
int useSimWt_
Definition: HGCalSD.h:56
std::unique_ptr< HGCGuardRing > guardRing_
Definition: HGCalSD.h:45
double sensorSizeOffset(bool reco) const
int levelT1_
Definition: HGCalSD.h:53
std::unique_ptr< HGCMouseBite > mouseBite_
Definition: HGCalSD.h:47
double eminHit_
Definition: HGCalSD.h:51
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:860
int levelTop(int ind=0) const
static int32_t layerFrontBack(int32_t layerOrient)
Definition: HGCalTypes.h:125
double waferSize(bool reco) const
int layerType(int lay) const
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCalSD.cc:244