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