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