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