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 
59  if (storeAllG4Hits_) {
60  setUseMap(false);
62  }
63 
64  //this is defined in the hgcsens.xml
65  G4String myName = name;
67  nameX_ = "HGCal";
68  if (myName.find("HitsEE") != std::string::npos) {
70  nameX_ = "HGCalEESensitive";
71  } else if (myName.find("HitsHEfront") != std::string::npos) {
73  nameX_ = "HGCalHESiliconSensitive";
74  }
75 
76 #ifdef EDM_ML_DEBUG
77  edm::LogVerbatim("HGCSim") << "**************************************************"
78  << "\n"
79  << "* *"
80  << "\n"
81  << "* Constructing a HGCalSD with name " << name << "\n"
82  << "* *"
83  << "\n"
84  << "**************************************************";
85 #endif
86  edm::LogVerbatim("HGCSim") << "HGCalSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
87  << " detector " << mydet_;
88  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
89  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
90  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
91  edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_ << " cuts along " << angles_.size()
92  << " axes: " << angles_[0] << ", " << angles_[1];
93 }
94 
95 double HGCalSD::getEnergyDeposit(const G4Step* aStep) {
96  double r = aStep->GetPreStepPoint()->GetPosition().perp();
97  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
98 #ifdef EDM_ML_DEBUG
99  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
100  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
101  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
102  edm::LogVerbatim("HGCSim") << "HGCalSD: Hit from standard path from " << lv->GetName() << " for Track "
103  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
104  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
105 #endif
106  // Apply fiducial cut
107  if (r < z * slopeMin_) {
108 #ifdef EDM_ML_DEBUG
109  edm::LogVerbatim("HGCSim") << "HGCalSD: Fiducial Volume cut";
110 #endif
111  return 0.0;
112  }
113 
114  double wt1 = getResponseWt(aStep->GetTrack());
115  double wt2 = aStep->GetTrack()->GetWeight();
116  double destep = weight_ * wt1 * (aStep->GetTotalEnergyDeposit());
117  if (wt2 > 0)
118  destep *= wt2;
119 #ifdef EDM_ML_DEBUG
120  edm::LogVerbatim("HGCSim") << "HGCalSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << " Total weight "
121  << weight_ * wt1 * wt2 << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
122 #endif
123  return destep;
124 }
125 
126 uint32_t HGCalSD::setDetUnitId(const G4Step* aStep) {
127  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
128  const G4VTouchable* touch = preStepPoint->GetTouchable();
129 
130 #ifdef EDM_ML_DEBUG
131  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_;
132  printDetectorLevels(touch);
133 #endif
134  //determine the exact position in global coordinates in the mass geometry
135  G4ThreeVector hitPoint = preStepPoint->GetPosition();
136  float globalZ = touch->GetTranslation(0).z();
137  int iz(globalZ > 0 ? 1 : -1);
138 
139  int layer(0), module(-1), cell(-1);
141  if (useSimWt_ > 0) {
142  layer = touch->GetReplicaNumber(2);
143  module = touch->GetReplicaNumber(1);
144  } else if (touch->GetHistoryDepth() > levelT2_) {
145  layer = touch->GetReplicaNumber(4);
146  module = touch->GetReplicaNumber(3);
147  cell = touch->GetReplicaNumber(1);
148  } else {
149  layer = touch->GetReplicaNumber(3);
150  module = touch->GetReplicaNumber(2);
151  }
152 #ifdef EDM_ML_DEBUG
153  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_ << ":"
154  << useSimWt_ << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell "
155  << layer << ":" << module << ":" << cell;
156  printDetectorLevels(touch);
157 #endif
158  } else if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
159  layer = touch->GetReplicaNumber(0);
160 #ifdef EDM_ML_DEBUG
161  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
162  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
163  << module << ":" << cell;
164 #endif
165  } else {
166  layer = touch->GetReplicaNumber(3);
167  module = touch->GetReplicaNumber(2);
168  cell = touch->GetReplicaNumber(1);
169 #ifdef EDM_ML_DEBUG
170  edm::LogVerbatim("HGCSim") << "DepthsInside: " << touch->GetHistoryDepth() << " name "
171  << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":" << module
172  << ":" << cell;
173 #endif
174  }
175 #ifdef EDM_ML_DEBUG
176  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
177  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
178  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
179  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
180  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
181  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
182  << touch->GetReplicaNumber(4) << " "
183  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
184  << mat->GetName() << ":" << mat->GetRadlen();
185 #endif
186  // The following statement should be examined later before elimination
187  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
188  return 0;
189 
190  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
191  if (rejectMB_ && id != 0) {
192  auto uv = HGCSiliconDetId(id).waferUV();
193 #ifdef EDM_ML_DEBUG
194  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
195 #endif
196  if (mouseBite_->exclude(hitPoint, iz, uv.first, uv.second)) {
197  id = 0;
198 #ifdef EDM_ML_DEBUG
199  edm::LogVerbatim("HGCSim") << "Rejected by mousebite cutoff *****";
200 #endif
201  }
202  }
203 #ifdef EDM_ML_DEBUG
204  if (id != 0)
205  edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
206 #endif
207  return id;
208 }
209 
210 void HGCalSD::update(const BeginOfJob* job) {
211  if (hgcons_ != nullptr) {
214  levelT1_ = hgcons_->levelTop(0);
215  levelT2_ = hgcons_->levelTop(1);
217  double waferSize = hgcons_->waferSize(false);
218  double mouseBite = hgcons_->mouseBite(false);
219  mouseBiteCut_ = waferSize * tan30deg_ - mouseBite;
220 #ifdef EDM_ML_DEBUG
221  edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
222  << " top Level " << levelT1_ << ":" << levelT2_ << " useSimWt " << useSimWt_ << " wafer "
223  << waferSize << ":" << mouseBite;
224 #endif
225 
226  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_, missingFile_);
227  if (rejectMB_)
228  mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
229  } else {
230  throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
231  }
232 }
233 
235 
236 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
237  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
238 }
239 
240 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
241  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
242  if (cornerMinMask_ > 2) {
243  if (hgcons_->maskCell(DetId(id), cornerMinMask_)) {
244  id = 0;
245  ignoreRejection();
246  }
247  }
248  if (hgcons_->waferHexagon8File() || (id == 0))
249  ignoreRejection();
250  return id;
251 }
252 
253 bool HGCalSD::isItinFidVolume(const G4ThreeVector& pos) {
254  if (fiducialCut_) {
255  return (hgcons_->distFromEdgeHex(pos.x(), pos.y(), pos.z()) > distanceFromEdge_);
256  } else {
257  return true;
258  }
259 }
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
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
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:234
double getEnergyDeposit(const G4Step *) override
Definition: HGCalSD.cc:95
bool rejectMB_
Definition: HGCalSD.h:52
DetId::Detector mydet_
Definition: HGCalSD.h:45
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCalSD.cc:126
std::string missingFile_
Definition: HGCalSD.h:56
bool storeAllG4Hits_
Definition: HGCalSD.h:51
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
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:253
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
bool waferHexagon8File() const
double weight_
Definition: HGCalSD.h:49
bool fiducialCut_
Definition: HGCalSD.h:52
void printDetectorLevels(const G4VTouchable *) const
Definition: CaloSD.cc:1081
HLT enums.
bool filterHit(CaloG4Hit *, double) override
Definition: HGCalSD.cc:236
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:858
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:210