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