CMS 3D CMS Logo

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