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