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  guardRing_.reset(nullptr);
47  mouseBite_.reset(nullptr);
48 
49  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCSD");
50  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
51  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
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), moduleLev(-1), cell(-1);
142  if (useSimWt_ > 0) {
143  layer = touch->GetReplicaNumber(2);
144  moduleLev = 1;
145  } else if (touch->GetHistoryDepth() > levelT2_) {
146  layer = touch->GetReplicaNumber(4);
147  cell = touch->GetReplicaNumber(1);
148  moduleLev = 3;
149  } else {
150  layer = touch->GetReplicaNumber(3);
151  moduleLev = 2;
152  }
153  int module = touch->GetReplicaNumber(moduleLev);
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 << ":" << moduleLev << ":" << module << ":" << cell;
158  printDetectorLevels(touch);
159  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
160  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
161  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
162  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
163  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
164  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
165  << touch->GetReplicaNumber(4) << " "
166  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
167  << mat->GetName() << ":" << mat->GetRadlen();
168 #endif
169  // The following statement should be examined later before elimination
170  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
171  return 0;
172 
173  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
174  if ((rejectMB_ || fiducialCut_) && id != 0) {
175  auto uv = HGCSiliconDetId(id).waferUV();
176 #ifdef EDM_ML_DEBUG
177  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
178 #endif
179  G4ThreeVector local = (touch->GetHistory()->GetTransform(moduleLev).TransformPoint(hitPoint));
180  if (fiducialCut_) {
181  int layertype = hgcons_->layerType(layer);
182  int frontBack = HGCalTypes::layerFrontBack(layertype);
183  if (guardRing_->exclude(local, iz, frontBack, layer, uv.first, uv.second)) {
184  id = 0;
185 #ifdef EDM_ML_DEBUG
186  edm::LogVerbatim("HGCSim") << "Rejected by GuardRing cutoff *****";
187 #endif
188  }
189  }
190  if ((rejectMB_) && (mouseBite_->exclude(local, iz, layer, uv.first, uv.second))) {
191  id = 0;
192 #ifdef EDM_ML_DEBUG
193  edm::LogVerbatim("HGCSim") << "Rejected by MouseBite cutoff *****";
194 #endif
195  }
196  }
197 #ifdef EDM_ML_DEBUG
198  if (id != 0)
199  edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
200 #endif
201  if ((id != 0) && checkID_) {
202  HGCSiliconDetId hid1(id);
203  bool cshift = (hgcons_->cassetteShiftSilicon(hid1.layer(), hid1.waferU(), hid1.waferV()));
204  std::string_view pid = (cshift ? "HGCSim" : "HGCalSim");
205  bool debug = (verbose_ > 0) ? true : false;
206  auto xy = hgcons_->locateCell(hid1, debug);
207  double xx = (hid1.zside() > 0) ? xy.first : -xy.first;
208  double dx = xx - (hitPoint.x() / CLHEP::cm);
209  double dy = xy.second - (hitPoint.y() / CLHEP::cm);
210  double diff = std::sqrt(dx * dx + dy * dy);
211  constexpr double tol = 2.0;
212  bool valid1 = hgcons_->isValidHex8(hid1.layer(), hid1.waferU(), hid1.waferV(), hid1.cellU(), hid1.cellV(), true);
213  if ((diff > tol) || (!valid1))
214  pid = "HGCalError";
215  auto partn = hgcons_->waferTypeRotation(hid1.layer(), hid1.waferU(), hid1.waferV(), false, false);
216  edm::LogVerbatim(pid) << "CheckID " << HGCSiliconDetId(id) << " input position: (" << hitPoint.x() / CLHEP::cm
217  << ", " << hitPoint.y() / CLHEP::cm << "); position from ID (" << xx << ", " << xy.second
218  << ") distance " << diff << " Valid " << valid1 << " Wafer type|rotation " << partn.first
219  << ":" << partn.second << " CassetteShift " << cshift;
220  }
221  return id;
222 }
223 
224 void HGCalSD::update(const BeginOfJob* job) {
225  if (hgcons_ != nullptr) {
228  levelT1_ = hgcons_->levelTop(0);
229  levelT2_ = hgcons_->levelTop(1);
231  int useOffset = hgcons_->getParameter()->useOffset_;
232  double waferSize = hgcons_->waferSize(false);
233  double mouseBite = hgcons_->mouseBite(false);
234  mouseBiteCut_ = waferSize * tan30deg_ - mouseBite;
235  if (useOffset > 0) {
236  rejectMB_ = true;
237  fiducialCut_ = true;
238  }
239 #ifdef EDM_ML_DEBUG
240  edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
241  << " top Level " << levelT1_ << ":" << levelT2_ << " useSimWt " << useSimWt_ << " wafer "
242  << waferSize << ":" << mouseBite << " useOffset " << useOffset;
243 #endif
244 
245  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_, missingFile_);
246  if (rejectMB_)
247  mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
248  if (fiducialCut_)
249  guardRing_ = std::make_unique<HGCGuardRing>(*hgcons_);
250  } else {
251  throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
252  }
253 }
254 
256 
257 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
258  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
259 }
260 
261 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
262  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
263  if (cornerMinMask_ > 2) {
264  if (hgcons_->maskCell(DetId(id), cornerMinMask_)) {
265  id = 0;
266  ignoreRejection();
267  }
268  }
269  if (hgcons_->waferHexagon8File() || (id == 0))
270  ignoreRejection();
271  return id;
272 }
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:28
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:255
double getEnergyDeposit(const G4Step *) override
Definition: HGCalSD.cc:97
bool rejectMB_
Definition: HGCalSD.h:54
DetId::Detector mydet_
Definition: HGCalSD.h:47
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:58
bool storeAllG4Hits_
Definition: HGCalSD.h:53
bool cassetteShiftSilicon(int layer, int waferU, int waferV) const
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
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:257
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:224