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 
16 #include "G4LogicalVolumeStore.hh"
17 #include "G4LogicalVolume.hh"
18 #include "G4Step.hh"
19 #include "G4Track.hh"
20 #include "G4ParticleTable.hh"
21 #include "G4VProcess.hh"
22 #include "G4Trap.hh"
23 
24 #include <iostream>
25 #include <fstream>
26 #include <iomanip>
27 
28 //#define EDM_ML_DEBUG
29 
31  const DDCompactView& cpv,
32  const SensitiveDetectorCatalog& clg,
33  edm::ParameterSet const& p,
34  const SimTrackManager* manager)
35  : CaloSD(name,
36  cpv,
37  clg,
38  p,
39  manager,
40  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
41  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
42  hgcons_(nullptr),
43  slopeMin_(0),
44  levelT1_(99),
45  levelT2_(99),
46  tan30deg_(std::tan(30.0 * CLHEP::deg)) {
47  numberingScheme_.reset(nullptr);
48  mouseBite_.reset(nullptr);
49 
51  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
52  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
53  distanceFromEdge_ = m_HGC.getParameter<double>("DistanceFromEdge");
54  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
55  rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
56  waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
57  cornerMinMask_ = m_HGC.getParameter<int>("CornerMinMask");
58  angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
59 
60  if (storeAllG4Hits_) {
61  setUseMap(false);
63  }
64 
65  //this is defined in the hgcsens.xml
66  G4String myName = name;
68  nameX_ = "HGCal";
69  if (myName.find("HitsEE") != std::string::npos) {
71  nameX_ = "HGCalEESensitive";
72  } else if (myName.find("HitsHEfront") != std::string::npos) {
74  nameX_ = "HGCalHESiliconSensitive";
75  }
76 
77 #ifdef EDM_ML_DEBUG
78  edm::LogVerbatim("HGCSim") << "**************************************************"
79  << "\n"
80  << "* *"
81  << "\n"
82  << "* Constructing a HGCalSD with name " << name << "\n"
83  << "* *"
84  << "\n"
85  << "**************************************************";
86 #endif
87  edm::LogVerbatim("HGCSim") << "HGCalSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
88  << " detector " << mydet_;
89  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
90  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
91  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
92  edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_ << " cuts along " << angles_.size()
93  << " axes: " << angles_[0] << ", " << angles_[1];
94 }
95 
96 double HGCalSD::getEnergyDeposit(const G4Step* aStep) {
97  double r = aStep->GetPreStepPoint()->GetPosition().perp();
98  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
99 #ifdef EDM_ML_DEBUG
100  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
101  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
102  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
103  edm::LogVerbatim("HGCSim") << "HGCalSD: Hit from standard path from " << lv->GetName() << " for Track "
104  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
105  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
106 #endif
107  // Apply fiducial cut
108  if (r < z * slopeMin_) {
109 #ifdef EDM_ML_DEBUG
110  edm::LogVerbatim("HGCSim") << "HGCalSD: Fiducial Volume cut";
111 #endif
112  return 0.0;
113  }
114 
115  double wt1 = getResponseWt(aStep->GetTrack());
116  double wt2 = aStep->GetTrack()->GetWeight();
117  double destep = weight_ * wt1 * (aStep->GetTotalEnergyDeposit());
118  if (wt2 > 0)
119  destep *= wt2;
120 #ifdef EDM_ML_DEBUG
121  edm::LogVerbatim("HGCSim") << "HGCalSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << " Total weight "
122  << weight_ * wt1 * wt2 << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
123 #endif
124  return destep;
125 }
126 
127 uint32_t HGCalSD::setDetUnitId(const G4Step* aStep) {
128  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
129  const G4VTouchable* touch = preStepPoint->GetTouchable();
130 
131  //determine the exact position in global coordinates in the mass geometry
132  G4ThreeVector hitPoint = preStepPoint->GetPosition();
133  float globalZ = touch->GetTranslation(0).z();
134  int iz(globalZ > 0 ? 1 : -1);
135 
136  int layer, module, cell;
137  if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
138  layer = touch->GetReplicaNumber(0);
139  module = -1;
140  cell = -1;
141 #ifdef EDM_ML_DEBUG
142  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
143  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
144  << module << ":" << cell;
145 #endif
146  } else {
147  layer = touch->GetReplicaNumber(3);
148  module = touch->GetReplicaNumber(2);
149  cell = touch->GetReplicaNumber(1);
150 #ifdef EDM_ML_DEBUG
151  edm::LogVerbatim("HGCSim") << "DepthsInside: " << touch->GetHistoryDepth() << " name "
152  << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":" << module
153  << ":" << cell;
154 #endif
155  }
156 #ifdef EDM_ML_DEBUG
157  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
158  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
159  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
160  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
161  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
162  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
163  << touch->GetReplicaNumber(4) << " "
164  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
165  << mat->GetName() << ":" << mat->GetRadlen();
166 #endif
167  // The following statement should be examined later before elimination
168  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
169  return 0;
170 
171  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
172  if (rejectMB_ && id != 0) {
173  auto uv = HGCSiliconDetId(id).waferUV();
174 #ifdef EDM_ML_DEBUG
175  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCSiliconDetId(id);
176 #endif
177  if (mouseBite_->exclude(hitPoint, iz, uv.first, uv.second)) {
178  id = 0;
179 #ifdef EDM_ML_DEBUG
180  edm::LogVerbatim("HGCSim") << "Rejected by mousebite cutoff *****";
181 #endif
182  }
183  }
184 #ifdef EDM_ML_DEBUG
185  if (id != 0)
186  edm::LogVerbatim("HGCSim") << HGCSiliconDetId(id);
187 #endif
188  return id;
189 }
190 
191 void HGCalSD::update(const BeginOfJob* job) {
192  const edm::EventSetup* es = (*job)();
194  es->get<IdealGeometryRecord>().get(nameX_, hdc);
195  if (hdc.isValid()) {
196  hgcons_ = hdc.product();
199  levelT1_ = hgcons_->levelTop(0);
200  levelT2_ = hgcons_->levelTop(1);
201  double waferSize = hgcons_->waferSize(false);
202  double mouseBite = hgcons_->mouseBite(false);
203  mouseBiteCut_ = waferSize * tan30deg_ - mouseBite;
204 #ifdef EDM_ML_DEBUG
205  edm::LogVerbatim("HGCSim") << "HGCalSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
206  << " top Level " << levelT1_ << ":" << levelT2_ << " wafer " << waferSize << ":"
207  << mouseBite;
208 #endif
209 
211  if (rejectMB_)
213  } else {
214  throw cms::Exception("Unknown", "HGCalSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
215  }
216 }
217 
219 
220 bool HGCalSD::filterHit(CaloG4Hit* aHit, double time) {
221  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
222 }
223 
224 uint32_t HGCalSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
225  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
226  if (cornerMinMask_ > 2) {
228  id = 0;
229  }
230  return id;
231 }
232 
233 bool HGCalSD::isItinFidVolume(const G4ThreeVector& pos) {
234  if (fiducialCut_) {
235  return (hgcons_->distFromEdgeHex(pos.x(), pos.y(), pos.z()) > distanceFromEdge_);
236  } else {
237  return true;
238  }
239 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
double distanceFromEdge_
Definition: HGCalSD.h:49
int levelT2_
Definition: HGCalSD.h:51
void setNumberCheckedHits(int val)
Definition: CaloSD.h:112
double slopeMin_
Definition: HGCalSD.h:49
Definition: CaloSD.h:38
#define nullptr
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCalSD.h:44
std::string nameX_
Definition: HGCalSD.h:47
double mouseBiteCut_
Definition: HGCalSD.h:50
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
HGCalSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCalSD.cc:30
void setUseMap(bool val)
Definition: CaloSD.h:101
const double MeV
void initRun() override
Definition: HGCalSD.cc:218
double getEnergyDeposit(const G4Step *) override
Definition: HGCalSD.cc:96
bool rejectMB_
Definition: HGCalSD.h:53
DetId::Detector mydet_
Definition: HGCalSD.h:46
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCalSD.cc:127
double mouseBite(bool reco) const
bool storeAllG4Hits_
Definition: HGCalSD.h:52
std::vector< double > angles_
Definition: HGCalSD.h:55
const HGCalDDDConstants * hgcons_
Definition: HGCalSD.h:43
double minSlope() const
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
HGCalGeometryMode::GeometryMode geomMode() const
double waferSize(bool reco) const
bool waferRot_
Definition: HGCalSD.h:53
std::pair< int, int > waferUV() const
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCalSD.cc:233
double tmaxHit
Definition: CaloSD.h:133
Definition: DetId.h:18
int cornerMinMask_
Definition: HGCalSD.h:51
bool maskCell(const DetId &id, int corners) const
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCalSD.h:48
double weight_
Definition: HGCalSD.h:50
bool fiducialCut_
Definition: HGCalSD.h:53
HLT enums.
bool filterHit(CaloG4Hit *, double) override
Definition: HGCalSD.cc:220
T get() const
Definition: EventSetup.h:71
int levelT1_
Definition: HGCalSD.h:51
std::unique_ptr< HGCMouseBite > mouseBite_
Definition: HGCalSD.h:45
int levelTop(int ind=0) const
double eminHit_
Definition: HGCalSD.h:49
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:623
Definition: vlib.h:208
double distFromEdgeHex(double x, double y, double z) const
double getEnergyDeposit() const
Definition: CaloG4Hit.h:77
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCalSD.cc:191