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