CMS 3D CMS Logo

HGCScintSD.cc
Go to the documentation of this file.
1 // File: HGCScintSD.cc
3 // Description: Sensitive Detector class for the Scintillator part of
4 // High Granularity Calorimeter
6 
17 #include "G4LogicalVolumeStore.hh"
18 #include "G4LogicalVolume.hh"
19 #include "G4Step.hh"
20 #include "G4Track.hh"
21 #include "G4ParticleTable.hh"
22 #include "G4VProcess.hh"
23 #include "G4Trap.hh"
24 
25 #include <fstream>
26 #include <iomanip>
27 #include <iostream>
28 #include <memory>
29 
30 //#define EDM_ML_DEBUG
31 
33  const HGCalDDDConstants* hgc,
34  const SensitiveDetectorCatalog& clg,
35  edm::ParameterSet const& p,
36  const SimTrackManager* manager)
37  : CaloSD(name,
38  clg,
39  p,
40  manager,
41  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
42  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
43  hgcons_(hgc),
44  slopeMin_(0),
45  levelT1_(99),
46  levelT2_(99),
47  firstLayer_(0) {
48  numberingScheme_.reset(nullptr);
49 
50  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCScintSD");
51  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
52  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
53  distanceFromEdge_ = m_HGC.getParameter<double>("DistanceFromEdge");
54  useBirk_ = m_HGC.getParameter<bool>("UseBirkLaw");
55  birk1_ = m_HGC.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
56  birk2_ = m_HGC.getParameter<double>("BirkC2");
57  birk3_ = m_HGC.getParameter<double>("BirkC3");
58  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
59  fileName_ = m_HGC.getUntrackedParameter<std::string>("TileFileName");
60 
61  if (storeAllG4Hits_) {
62  setUseMap(false);
64  }
65 
66  //this is defined in the hgcsens.xml
67  G4String myName = name;
69  nameX_ = "HGCal";
70  if (myName.find("HitsHEback") != std::string::npos) {
72  nameX_ = "HGCalHEScintillatorSensitive";
73  }
74 
75 #ifdef EDM_ML_DEBUG
76  edm::LogVerbatim("HGCSim") << "**************************************************"
77  << "\n"
78  << "* *"
79  << "\n"
80  << "* Constructing a HGCScintSD with name " << name << "\n"
81  << "* *"
82  << "\n"
83  << "**************************************************";
84 #endif
85  edm::LogVerbatim("HGCSim") << "HGCScintSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
86  << " detector " << mydet_ << " File " << fileName_;
87  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
88  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
89  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
90  edm::LogVerbatim("HGCSim") << "Use of Birks law is set to " << useBirk_
91  << " with three constants kB = " << birk1_ << ", C1 = " << birk2_ << ", C2 = " << birk3_;
92 
93  if (!fileName_.empty()) {
94  edm::FileInPath filetmp("SimG4CMS/Calo/data/" + fileName_);
95  std::string fileName = filetmp.fullPath();
96  std::ifstream fInput(fileName.c_str());
97  if (!fInput.good()) {
98  edm::LogVerbatim("HGCSim") << "Cannot open file " << fileName;
99  } else {
100  char buffer[80];
101  while (fInput.getline(buffer, 80)) {
102  std::vector<std::string> items = CaloSimUtils::splitString(std::string(buffer));
103  if (items.size() > 2) {
104  int layer = std::atoi(items[0].c_str());
105  int ring = std::atoi(items[1].c_str());
106  int phi = std::atoi(items[2].c_str());
107  tiles_.emplace_back(HGCalTileIndex::tileIndex(layer, ring, phi));
108  }
109  }
110  edm::LogVerbatim("HGCSim") << "Reads in " << tiles_.size() << " tile information from " << fileName_;
111  fInput.close();
112  }
113  }
114 }
115 
116 double HGCScintSD::getEnergyDeposit(const G4Step* aStep) {
117  double r = aStep->GetPreStepPoint()->GetPosition().perp();
118  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
119 #ifdef EDM_ML_DEBUG
120  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
121  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
122  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
123  edm::LogVerbatim("HGCSim") << "HGCScintSD: Hit from standard path from " << lv->GetName() << " for Track "
124  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
125  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
126 #endif
127  // Apply fiducial cut
128  if (r < z * slopeMin_) {
129 #ifdef EDM_ML_DEBUG
130  edm::LogVerbatim("HGCSim") << "HGCScintSD: Fiducial Volume cut";
131 #endif
132  return 0.0;
133  }
134 
135  double wt1 = getResponseWt(aStep->GetTrack());
136  double wt2 = aStep->GetTrack()->GetWeight();
137  double wt3 = (useBirk_ ? getAttenuation(aStep, birk1_, birk2_, birk3_) : 1.0);
138  double destep = weight_ * wt1 * wt3 * (aStep->GetTotalEnergyDeposit());
139  if (wt2 > 0)
140  destep *= wt2;
141 #ifdef EDM_ML_DEBUG
142  edm::LogVerbatim("HGCalSim") << "HGCScintSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << ":" << wt3
143  << " Total weight " << weight_ * wt1 * wt2 * wt3
144  << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
145 #endif
146  return destep;
147 }
148 
149 uint32_t HGCScintSD::setDetUnitId(const G4Step* aStep) {
150  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
151  const G4VTouchable* touch = preStepPoint->GetTouchable();
152 
153 #ifdef EDM_ML_DEBUG
154  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_;
155  printDetectorLevels(touch);
156 #endif
157  //determine the exact position in global coordinates in the mass geometry
158  G4ThreeVector hitPoint = preStepPoint->GetPosition();
159  float globalZ = touch->GetTranslation(0).z();
160  int iz(globalZ > 0 ? 1 : -1);
161 
162  int layer(0), module(-1), cell(-1);
164  layer = touch->GetReplicaNumber(1);
165  } else if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
166  layer = touch->GetReplicaNumber(0);
167 #ifdef EDM_ML_DEBUG
168  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
169  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
170  << module << ":" << cell;
171 #endif
172  } else {
173  layer = touch->GetReplicaNumber(3);
174  module = touch->GetReplicaNumber(2);
175  cell = touch->GetReplicaNumber(1);
176 #ifdef EDM_ML_DEBUG
177  edm::LogVerbatim("HGCSim") << "DepthsInside: " << touch->GetHistoryDepth() << " name "
178  << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":" << module
179  << ":" << cell;
180 #endif
181  }
182 #ifdef EDM_ML_DEBUG
183  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
184  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
185  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
186  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
187  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
188  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
189  << touch->GetReplicaNumber(4) << " "
190  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
191  << mat->GetName() << ":" << mat->GetRadlen();
192 #endif
193  // The following statement should be examined later before elimination
194  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
195  return 0;
196 
197  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
198  bool debug(false);
199  if (!tiles_.empty()) {
200  HGCScintillatorDetId hid(id);
201  int indx = HGCalTileIndex::tileIndex(firstLayer_ + hid.layer(), hid.ring(), hid.iphi());
202  if (std::find(tiles_.begin(), tiles_.end(), indx) != tiles_.end())
203  debug = true;
204  }
205  if (debug)
206  edm::LogVerbatim("HGCSim") << "Layer:module:cell:iz " << layer << ":" << module << ":" << cell << ":" << iz
207  << " Point (" << hitPoint.x() << ", " << hitPoint.y() << ", " << hitPoint.z() << ") "
208  << HGCScintillatorDetId(id);
209 
210  if (!isItinFidVolume(hitPoint)) {
211 #ifdef EDM_ML_DEBUG
212  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCScintillatorDetId(id)
213  << " is rejected by fiducilal volume cut";
214 #endif
215  id = 0;
216  }
217  return id;
218 }
219 
220 void HGCScintSD::update(const BeginOfJob* job) {
221  if (hgcons_ != nullptr) {
224  levelT1_ = hgcons_->levelTop(0);
225  levelT2_ = hgcons_->levelTop(1);
226  firstLayer_ = hgcons_->firstLayer() - 1;
227 #ifdef EDM_ML_DEBUG
228  edm::LogVerbatim("HGCSim") << "HGCScintSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
229  << " top Level " << levelT1_ << ":" << levelT2_ << " FirstLayer " << firstLayer_;
230 #endif
231 
232  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_, fileName_);
233  } else {
234  throw cms::Exception("Unknown", "HGCScintSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
235  }
236 }
237 
239 
240 bool HGCScintSD::filterHit(CaloG4Hit* aHit, double time) {
241  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
242 }
243 
244 uint32_t HGCScintSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
245  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
246  return id;
247 }
248 
249 bool HGCScintSD::isItinFidVolume(const G4ThreeVector& pos) {
250  if (fiducialCut_) {
251  return (hgcons_->distFromEdgeTrap(pos.x(), pos.y(), pos.z()) > distanceFromEdge_);
252  } else {
253  return true;
254  }
255 }
double birk2_
Definition: HGCScintSD.h:50
Log< level::Info, true > LogVerbatim
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCScintSD.h:45
int levelT2_
Definition: HGCScintSD.h:47
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void setNumberCheckedHits(int val)
Definition: CaloSD.h:122
double birk1_
Definition: HGCScintSD.h:50
HGCScintSD(const std::string &, const HGCalDDDConstants *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCScintSD.cc:32
Definition: CaloSD.h:40
std::string fullPath() const
Definition: FileInPath.cc:161
double distFromEdgeTrap(double x, double y, double z) const
double distanceFromEdge_
Definition: HGCScintSD.h:46
std::vector< int > tiles_
Definition: HGCScintSD.h:52
int firstLayer() const
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
HGCalGeometryMode::GeometryMode geomMode() const
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCScintSD.cc:149
int levelT1_
Definition: HGCScintSD.h:47
void setUseMap(bool val)
Definition: CaloSD.h:111
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
bool fiducialCut_
Definition: HGCScintSD.h:48
bool useBirk_
Definition: HGCScintSD.h:49
constexpr std::array< uint8_t, layerIndexSize > layer
const HGCalDDDConstants * hgcons_
Definition: HGCScintSD.h:41
T getUntrackedParameter(std::string const &, T const &) const
int32_t tileIndex(int32_t layer, int32_t ring, int32_t phi)
int iphi() const
get the phi index
std::vector< std::string > splitString(const std::string &)
Definition: CaloSimUtils.cc:3
double weight_
Definition: HGCScintSD.h:50
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCScintSD.cc:249
int layer() const
get the layer #
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void initRun() override
Definition: HGCScintSD.cc:238
double tmaxHit
Definition: CaloSD.h:144
int ring() const
get the eta index
double getEnergyDeposit(const G4Step *) override
Definition: HGCScintSD.cc:116
#define debug
Definition: HDRShower.cc:19
double eminHit_
Definition: HGCScintSD.h:46
double minSlope() const
std::string fileName_
Definition: HGCScintSD.h:51
double slopeMin_
Definition: HGCScintSD.h:46
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:666
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCScintSD.h:42
int firstLayer_
Definition: HGCScintSD.h:47
void printDetectorLevels(const G4VTouchable *) const
Definition: CaloSD.cc:1081
std::string nameX_
Definition: HGCScintSD.h:44
HLT enums.
DetId::Detector mydet_
Definition: HGCScintSD.h:43
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:858
int levelTop(int ind=0) const
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCScintSD.cc:220
double birk3_
Definition: HGCScintSD.h:50
bool storeAllG4Hits_
Definition: HGCScintSD.h:48
bool filterHit(CaloG4Hit *, double) override
Definition: HGCScintSD.cc:240