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  angles_ = m_HFN.getUntrackedParameter<std::vector<double>>("WaferAngles");
49 
50  nameX_ = ((name.find("HFNoseHits")!=std::string::npos) ?
51  "HGCalHFNoseSensitive" : "HFNoseSensitive");
52 
53  if(storeAllG4Hits_) {
54  setUseMap(false);
56  }
57 
58 #ifdef EDM_ML_DEBUG
59  edm::LogVerbatim("HFNSim")<< "**************************************************"
60  << "\n"
61  << "* *"
62  << "\n"
63  << "* Constructing a HFNoseSD with name " << name << "\n"
64  << "* *"
65  << "\n"
66  << "**************************************************";
67 #endif
68  edm::LogVerbatim("HFNSim") << "HFNoseSD:: Threshold for storing hits: "
69  << eminHit_ << " for " << name;
70  edm::LogVerbatim("HFNSim") << "Flag for storing individual Geant4 Hits "
71  << storeAllG4Hits_;
72  edm::LogVerbatim("HFNSim") << "Fiducial volume cut with cut from eta/phi "
73  << "boundary " << fiducialCut_ << " at "
75  edm::LogVerbatim("HFNSim") << "Reject MosueBite Flag: " << rejectMB_
76  << " cuts along " << angles_.size() << " axes: "
77  << angles_[0] << ", " << angles_[1];
78 }
79 
80 double HFNoseSD::getEnergyDeposit(const G4Step* aStep) {
81 
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 =
88  aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
89  edm::LogVerbatim("HFNSim") << "HFNoseSD: Hit from standard path from "
90  << lv->GetName() << " for Track "
91  << aStep->GetTrack()->GetTrackID() << " ("
92  << parCode << ":" << parName << ") R = " << r
93  << " Z = " << z << " slope = " << r/z << ":"
94  << slopeMin_;
95 #endif
96  // Apply fiducial cut
97  if (r < z*slopeMin_) {
98 #ifdef EDM_ML_DEBUG
99  edm::LogVerbatim("HFNSim") << "HFNoseSD: Fiducial Volume cut";
100 #endif
101  return 0.0;
102  }
103 
104  double wt1 = getResponseWt(aStep->GetTrack());
105  double wt2 = aStep->GetTrack()->GetWeight();
106  double destep = weight_*wt1*(aStep->GetTotalEnergyDeposit());
107  if (wt2 > 0) destep *= wt2;
108 #ifdef EDM_ML_DEBUG
109  edm::LogVerbatim("HFNSim") << "HFNoseSD: weights= " << weight_ << ":"
110  << wt1 << ":" << wt2 << " Total weight "
111  << weight_*wt1*wt2 << " deStep: "
112  << aStep->GetTotalEnergyDeposit()
113  << ":" <<destep;
114 #endif
115  return destep;
116 }
117 
118 uint32_t HFNoseSD::setDetUnitId(const G4Step * aStep) {
119 
120  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
121  const G4VTouchable* touch = preStepPoint->GetTouchable();
122 
123  //determine the exact position in global coordinates in the mass geometry
124  G4ThreeVector hitPoint = preStepPoint->GetPosition();
125  float globalZ = touch->GetTranslation(0).z();
126  int iz( globalZ>0 ? 1 : -1);
127 
128  int layer, module, cell;
129  if ((touch->GetHistoryDepth() == levelT1_) ||
130  (touch->GetHistoryDepth() == levelT2_)) {
131  layer = touch->GetReplicaNumber(0);
132  module = -1;
133  cell = -1;
134 #ifdef EDM_ML_DEBUG
135  edm::LogVerbatim("HFNSim") << "DepthsTop: " << touch->GetHistoryDepth()
136  << ":" << levelT1_ << ":" << levelT2_
137  << " name " << touch->GetVolume(0)->GetName()
138  << " layer:module:cell " << layer << ":"
139  << module << ":" << cell;
140 #endif
141  } else {
142  layer = touch->GetReplicaNumber(3);
143  module = touch->GetReplicaNumber(2);
144  cell = touch->GetReplicaNumber(1);
145 #ifdef EDM_ML_DEBUG
146  edm::LogVerbatim("HFNSim") << "DepthsInside: " << touch->GetHistoryDepth()
147  << " name " << touch->GetVolume(0)->GetName()
148  << " layer:module:cell " << layer << ":"
149  << module << ":" << cell;
150 #endif
151  }
152 #ifdef EDM_ML_DEBUG
153  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
154  edm::LogVerbatim("HFNSim") << "Depths: " << touch->GetHistoryDepth()
155  << " name " << touch->GetVolume(0)->GetName()
156  << ":" << touch->GetReplicaNumber(0) << " "
157  << touch->GetVolume(1)->GetName()
158  << ":" << touch->GetReplicaNumber(1) << " "
159  << touch->GetVolume(2)->GetName()
160  << ":" << touch->GetReplicaNumber(2) << " "
161  << touch->GetVolume(3)->GetName()
162  << ":" << touch->GetReplicaNumber(3) << " "
163  << touch->GetVolume(4)->GetName()
164  << ":" << touch->GetReplicaNumber(4) << " "
165  << " layer:module:cell " << layer << ":"
166  << module << ":" << cell << " Material "
167  << mat->GetName() << ":" << mat->GetRadlen();
168 #endif
169  // The following statement should be examined later before elimination
170  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.) return 0;
171 
172  uint32_t id = setDetUnitId (layer, module, cell, iz, hitPoint);
173  if (rejectMB_ && id != 0) {
174  auto uv = HFNoseDetId(id).waferUV();
175 #ifdef EDM_ML_DEBUG
176  edm::LogVerbatim("HFNSim") << "ID " << std::hex << id << std::dec
177  << " " << HFNoseDetId(id);
178 #endif
179  if (mouseBite_->exclude(hitPoint, iz, uv.first, uv.second)) {
180  id = 0;
181 #ifdef EDM_ML_DEBUG
182  edm::LogVerbatim("HFNSim") << "Rejected by mousebite cutoff *****";
183 #endif
184  }
185  }
186  return id;
187 }
188 
189 void HFNoseSD::update(const BeginOfJob * job) {
190 
191  const edm::EventSetup* es = (*job)();
193  es->get<IdealGeometryRecord>().get(nameX_,hdc);
194  if (hdc.isValid()) {
195  hgcons_ = hdc.product();
198  levelT1_ = hgcons_->levelTop(0);
199  levelT2_ = hgcons_->levelTop(1);
200  double waferSize = hgcons_->waferSize(false);
201  double mouseBite = hgcons_->mouseBite(false);
202  mouseBiteCut_ = waferSize*tan30deg_ - mouseBite;
203 #ifdef EDM_ML_DEBUG
204  edm::LogVerbatim("HFNSim") << "HFNoseSD::Initialized with mode "
205  << geom_mode_ << " Slope cut " << slopeMin_
206  << " top Level " << levelT1_ << ":" << levelT2_
207  << " wafer " << waferSize << ":" << mouseBite;
208 #endif
209 
211  if (rejectMB_)
213  waferRot_));
214  } else {
215  throw cms::Exception("Unknown", "HFNoseSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
216  }
217 }
218 
220 }
221 
222 bool HFNoseSD::filterHit(CaloG4Hit* aHit, double time) {
223  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
224 }
225 
226 uint32_t HFNoseSD::setDetUnitId (int layer, int module, int cell, int iz,
227  G4ThreeVector &pos) {
228  uint32_t id = numberingScheme_ ?
229  numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
230  return id;
231 }
232 
233 bool HFNoseSD::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 eminHit_
Definition: HFNoseSD.h:50
std::unique_ptr< HFNoseNumberingScheme > numberingScheme_
Definition: HFNoseSD.h:46
void setNumberCheckedHits(int val)
Definition: CaloSD.h:111
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:118
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:233
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
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:189
double getEnergyDeposit(const G4Step *) override
Definition: HFNoseSD.cc:80
std::string nameX_
Definition: HFNoseSD.h:48
HLT enums.
void initRun() override
Definition: HFNoseSD.cc:219
T get() const
Definition: EventSetup.h:62
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:222
double slopeMin_
Definition: HFNoseSD.h:50