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  numberingScheme_.reset(nullptr);
49  guardRing_.reset(nullptr);
50  mouseBite_.reset(nullptr);
51 
52  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCSD");
53  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
54  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
55  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
56  rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
57  waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
58  cornerMinMask_ = m_HGC.getParameter<int>("CornerMinMask");
59  angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
60  missingFile_ = m_HGC.getUntrackedParameter<std::string>("MissingWaferFile");
61  checkID_ = m_HGC.getUntrackedParameter<bool>("CheckID");
62  verbose_ = m_HGC.getUntrackedParameter<int>("Verbosity");
63 
64  if (storeAllG4Hits_) {
65  setUseMap(false);
67  }
68 
69  //this is defined in the hgcsens.xml
70  G4String myName = name;
72  nameX_ = "HGCal";
73  if (myName.find("HitsEE") != std::string::npos) {
75  nameX_ = "HGCalEESensitive";
76  } else if (myName.find("HitsHEfront") != std::string::npos) {
78  nameX_ = "HGCalHESiliconSensitive";
79  }
80 
81 #ifdef EDM_ML_DEBUG
82  edm::LogVerbatim("HGCSim") << "**************************************************"
83  << "\n"
84  << "* *"
85  << "\n"
86  << "* Constructing a HGCalSD with name " << name << "\n"
87  << "* *"
88  << "\n"
89  << "**************************************************";
90 #endif
91  edm::LogVerbatim("HGCSim") << "HGCalSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
92  << " detector " << mydet_;
93  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
94  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
95  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
96  edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_ << " cuts along " << angles_.size()
97  << " axes: " << angles_[0] << ", " << angles_[1];
98 }
99 
100 double HGCalSD::getEnergyDeposit(const G4Step* aStep) {
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  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
106  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
107  edm::LogVerbatim("HGCSim") << "HGCalSD: Hit from standard path from " << lv->GetName() << " for Track "
108  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
109  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
110 #endif
111  // Apply fiducial cut
112  if (r < z * slopeMin_) {
113 #ifdef EDM_ML_DEBUG
114  edm::LogVerbatim("HGCSim") << "HGCalSD: Fiducial Volume cut";
115 #endif
116  return 0.0;
117  }
118 
119  double wt1 = getResponseWt(aStep->GetTrack());
120  double wt2 = aStep->GetTrack()->GetWeight();
121  double destep = weight_ * wt1 * (aStep->GetTotalEnergyDeposit());
122  if (wt2 > 0)
123  destep *= wt2;
124 #ifdef EDM_ML_DEBUG
125  edm::LogVerbatim("HGCSim") << "HGCalSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << " Total weight "
126  << weight_ * wt1 * wt2 << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
127 #endif
128  return destep;
129 }
130 
131 uint32_t HGCalSD::setDetUnitId(const G4Step* aStep) {
132  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
133  const G4VTouchable* touch = preStepPoint->GetTouchable();
134 
135 #ifdef EDM_ML_DEBUG
136  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_;
137  printDetectorLevels(touch);
138 #endif
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  int layer(0), moduleLev(-1), cell(-1);
145  if (useSimWt_ > 0) {
146  layer = touch->GetReplicaNumber(2);
147  moduleLev = 1;
148  } else if (touch->GetHistoryDepth() > levelT2_) {
149  layer = touch->GetReplicaNumber(4);
150  cell = touch->GetReplicaNumber(1);
151  moduleLev = 3;
152  } else {
153  layer = touch->GetReplicaNumber(3);
154  moduleLev = 2;
155  }
156  int module = touch->GetReplicaNumber(moduleLev);
157 #ifdef EDM_ML_DEBUG
158  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << ":"
159  << useSimWt_ << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell "
160  << layer << ":" << moduleLev << ":" << module << ":" << cell;
161  printDetectorLevels(touch);
162  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
163  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
164  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
165  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
166  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
167  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
168  << touch->GetReplicaNumber(4) << " "
169  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
170  << mat->GetName() << ":" << mat->GetRadlen();
171 #endif
172  // The following statement should be examined later before elimination
173  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
174  return 0;
175 
176  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
177  if ((rejectMB_ || fiducialCut_) && id != 0) {
178  auto uv = HGCSiliconDetId(id).waferUV();
179 #ifdef EDM_ML_DEBUG
180  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
181 #endif
182  G4ThreeVector local = (touch->GetHistory()->GetTransform(moduleLev).TransformPoint(hitPoint));
183  if (fiducialCut_) {
184  int layertype = hgcons_->layerType(layer);
185  int frontBack = HGCalTypes::layerFrontBack(layertype);
186  if (guardRing_->exclude(local, iz, frontBack, layer, uv.first, uv.second)) {
187  id = 0;
188 #ifdef EDM_ML_DEBUG
189  edm::LogVerbatim("HGCSim") << "Rejected by GuardRing cutoff *****";
190 #endif
191  }
192  }
193  if ((rejectMB_) && (mouseBite_->exclude(local, iz, layer, uv.first, uv.second))) {
194  id = 0;
195 #ifdef EDM_ML_DEBUG
196  edm::LogVerbatim("HGCSim") << "Rejected by MouseBite cutoff *****";
197 #endif
198  }
199  }
200 #ifdef EDM_ML_DEBUG
201  if (id != 0)
202  edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
203 #endif
204  if ((id != 0) && checkID_) {
205  HGCSiliconDetId hid1(id);
206  bool cshift = (hgcons_->cassetteShiftSilicon(hid1.zside(), hid1.layer(), hid1.waferU(), hid1.waferV()));
207  std::string_view pid = (cshift ? "HGCSim" : "HGCalSim");
208  bool debug = (verbose_ > 0) ? true : false;
209  auto xy = hgcons_->locateCell(hid1, debug);
210  double xx = (hid1.zside() > 0) ? xy.first : -xy.first;
211  double dx = xx - (hitPoint.x() / CLHEP::cm);
212  double dy = xy.second - (hitPoint.y() / CLHEP::cm);
213  double diff = (dx * dx + dy * dy);
214  constexpr double tol = 2.0 * 2.0;
215  bool valid1 = hgcons_->isValidHex8(hid1.layer(), hid1.waferU(), hid1.waferV(), hid1.cellU(), hid1.cellV(), true);
216  if ((diff > tol) || (!valid1))
217  pid = "HGCalError";
218  auto partn = hgcons_->waferTypeRotation(hid1.layer(), hid1.waferU(), hid1.waferV(), false, false);
219  int indx = HGCalWaferIndex::waferIndex(layer, hid1.waferU(), hid1.waferV());
220  edm::LogVerbatim(pid) << "CheckID " << HGCSiliconDetId(id) << " Layer:Module:Cell:ModuleLev " << layer << ":"
221  << module << ":" << cell << ":" << moduleLev << " SimWt:history " << useSimWt_ << ":"
222  << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << " input position: ("
223  << hitPoint.x() / CLHEP::cm << ", " << hitPoint.y() / CLHEP::cm << ":"
224  << convertRadToDeg(std::atan2(hitPoint.y(), hitPoint.x())) << "); position from ID (" << xx
225  << ", " << xy.second << ") distance " << dx << ":" << dy << ":" << std::sqrt(diff)
226  << " Valid " << valid1 << " Wafer type|rotation " << partn.first << ":" << partn.second
227  << " Part:Orient:Cassette " << std::get<1>(hgcons_->waferFileInfo(indx)) << ":"
228  << std::get<2>(hgcons_->waferFileInfo(indx)) << ":"
229  << std::get<3>(hgcons_->waferFileInfo(indx)) << " CassetteShift " << cshift;
230  if ((diff > tol) || (!valid1)) {
231  printDetectorLevels(touch);
232  hgcons_->locateCell(hid1, true);
233  }
234  }
235  return id;
236 }
237 
238 void HGCalSD::update(const BeginOfJob* job) {
239  if (hgcons_ != nullptr) {
242  levelT1_ = hgcons_->levelTop(0);
243  levelT2_ = hgcons_->levelTop(1);
245  int useOffset = hgcons_->getParameter()->useOffset_;
246  double waferSize = hgcons_->waferSize(false);
247  double mouseBite = hgcons_->mouseBite(false);
248  mouseBiteCut_ = waferSize * tan30deg_ - mouseBite;
249  if (useOffset > 0) {
250  rejectMB_ = true;
251  fiducialCut_ = true;
252  }
253 #ifdef EDM_ML_DEBUG
254  edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
255  << " top Level " << levelT1_ << ":" << levelT2_ << " useSimWt " << useSimWt_ << " wafer "
256  << waferSize << ":" << mouseBite << " useOffset " << useOffset;
257 #endif
258 
259  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_, missingFile_);
260  if (rejectMB_)
261  mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
262  if (fiducialCut_)
263  guardRing_ = std::make_unique<HGCGuardRing>(*hgcons_);
264  } else {
265  throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
266  }
267 }
268 
270 
271 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
272  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
273 }
274 
275 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
276  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
277  if (cornerMinMask_ > 2) {
278  if (hgcons_->maskCell(DetId(id), cornerMinMask_)) {
279  id = 0;
280  ignoreRejection();
281  }
282  }
283  if (hgcons_->waferHexagon8File() || (id == 0))
284  ignoreRejection();
285  return id;
286 }
bool maskCell(const DetId &id, int corners) const
Log< level::Info, true > LogVerbatim
int verbose_
Definition: HGCalSD.h:55
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double distanceFromEdge_
Definition: HGCalSD.h:50
int levelT2_
Definition: HGCalSD.h:52
void setNumberCheckedHits(int val)
Definition: CaloSD.h:122
double slopeMin_
Definition: HGCalSD.h:50
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:48
double mouseBiteCut_
Definition: HGCalSD.h:51
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
constexpr std::array< uint8_t, layerIndexSize< TrackerTraits > > layer
int layer() const
get the layer #
void initRun() override
Definition: HGCalSD.cc:269
double getEnergyDeposit(const G4Step *) override
Definition: HGCalSD.cc:100
bool rejectMB_
Definition: HGCalSD.h:54
DetId::Detector mydet_
Definition: HGCalSD.h:47
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCalSD.cc:131
std::tuple< int, int, int, int > waferFileInfo(unsigned int kk) const
T sqrt(T t)
Definition: SSEVec.h:19
std::string missingFile_
Definition: HGCalSD.h:58
bool storeAllG4Hits_
Definition: HGCalSD.h:53
std::vector< double > angles_
Definition: HGCalSD.h:57
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:56
double mouseBite(bool reco) const
bool waferRot_
Definition: HGCalSD.h:54
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:54
Definition: DetId.h:17
int cornerMinMask_
Definition: HGCalSD.h:52
#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:49
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:51
bool fiducialCut_
Definition: HGCalSD.h:54
void printDetectorLevels(const G4VTouchable *) const
Definition: CaloSD.cc:1078
HLT enums.
bool filterHit(CaloG4Hit *, double) override
Definition: HGCalSD.cc:271
bool cassetteShiftSilicon(int zside, int layer, int waferU, int waferV) const
int useSimWt_
Definition: HGCalSD.h:55
std::unique_ptr< HGCGuardRing > guardRing_
Definition: HGCalSD.h:45
int levelT1_
Definition: HGCalSD.h:52
std::unique_ptr< HGCMouseBite > mouseBite_
Definition: HGCalSD.h:46
double eminHit_
Definition: HGCalSD.h:50
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:855
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:238