CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
15 #include "G4LogicalVolumeStore.hh"
16 #include "G4LogicalVolume.hh"
17 #include "G4Step.hh"
18 #include "G4Track.hh"
19 #include "G4ParticleTable.hh"
20 #include "G4VProcess.hh"
21 #include "G4Trap.hh"
22 
23 #include <fstream>
24 #include <iomanip>
25 #include <iostream>
26 #include <memory>
27 
28 //#define EDM_ML_DEBUG
29 
31  const HGCalDDDConstants* hgc,
32  const SensitiveDetectorCatalog& clg,
33  edm::ParameterSet const& p,
34  const SimTrackManager* manager)
35  : CaloSD(name,
36  clg,
37  p,
38  manager,
39  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
40  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
41  hgcons_(hgc),
42  slopeMin_(0),
43  levelT1_(99),
44  levelT2_(99) {
45  numberingScheme_.reset(nullptr);
46 
47  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCScintSD");
48  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
49  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
50  distanceFromEdge_ = m_HGC.getParameter<double>("DistanceFromEdge");
51  useBirk_ = m_HGC.getParameter<bool>("UseBirkLaw");
52  birk1_ = m_HGC.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
53  birk2_ = m_HGC.getParameter<double>("BirkC2");
54  birk3_ = m_HGC.getParameter<double>("BirkC3");
55  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
56 
57  if (storeAllG4Hits_) {
58  setUseMap(false);
60  }
61 
62  //this is defined in the hgcsens.xml
63  G4String myName = name;
65  nameX_ = "HGCal";
66  if (myName.find("HitsHEback") != std::string::npos) {
68  nameX_ = "HGCalHEScintillatorSensitive";
69  }
70 
71 #ifdef EDM_ML_DEBUG
72  edm::LogVerbatim("HGCSim") << "**************************************************"
73  << "\n"
74  << "* *"
75  << "\n"
76  << "* Constructing a HGCScintSD with name " << name << "\n"
77  << "* *"
78  << "\n"
79  << "**************************************************";
80 #endif
81  edm::LogVerbatim("HGCSim") << "HGCScintSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
82  << " detector " << mydet_;
83  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
84  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
85  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
86  edm::LogVerbatim("HGCSim") << "Use of Birks law is set to " << useBirk_
87  << " with three constants kB = " << birk1_ << ", C1 = " << birk2_ << ", C2 = " << birk3_;
88 }
89 
90 double HGCScintSD::getEnergyDeposit(const G4Step* aStep) {
91  double r = aStep->GetPreStepPoint()->GetPosition().perp();
92  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
93 #ifdef EDM_ML_DEBUG
94  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
95  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
96  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
97  edm::LogVerbatim("HGCSim") << "HGCScintSD: Hit from standard path from " << lv->GetName() << " for Track "
98  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
99  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
100 #endif
101  // Apply fiducial cut
102  if (r < z * slopeMin_) {
103 #ifdef EDM_ML_DEBUG
104  edm::LogVerbatim("HGCSim") << "HGCScintSD: Fiducial Volume cut";
105 #endif
106  return 0.0;
107  }
108 
109  double wt1 = getResponseWt(aStep->GetTrack());
110  double wt2 = aStep->GetTrack()->GetWeight();
111  double wt3 = (useBirk_ ? getAttenuation(aStep, birk1_, birk2_, birk3_) : 1.0);
112  double destep = weight_ * wt1 * wt3 * (aStep->GetTotalEnergyDeposit());
113  if (wt2 > 0)
114  destep *= wt2;
115 #ifdef EDM_ML_DEBUG
116  edm::LogVerbatim("HGCalSim") << "HGCScintSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << ":" << wt3
117  << " Total weight " << weight_ * wt1 * wt2 * wt3
118  << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
119 #endif
120  return destep;
121 }
122 
123 uint32_t HGCScintSD::setDetUnitId(const G4Step* aStep) {
124  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
125  const G4VTouchable* touch = preStepPoint->GetTouchable();
126 
127 #ifdef EDM_ML_DEBUG
128  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_;
129  printDetectorLevels(touch);
130 #endif
131  //determine the exact position in global coordinates in the mass geometry
132  G4ThreeVector hitPoint = preStepPoint->GetPosition();
133  float globalZ = touch->GetTranslation(0).z();
134  int iz(globalZ > 0 ? 1 : -1);
135 
136  int layer(0), module(-1), cell(-1);
138  layer = touch->GetReplicaNumber(1);
139  } else if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
140  layer = touch->GetReplicaNumber(0);
141 #ifdef EDM_ML_DEBUG
142  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
143  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
144  << module << ":" << cell;
145 #endif
146  } else {
147  layer = touch->GetReplicaNumber(3);
148  module = touch->GetReplicaNumber(2);
149  cell = touch->GetReplicaNumber(1);
150 #ifdef EDM_ML_DEBUG
151  edm::LogVerbatim("HGCSim") << "DepthsInside: " << touch->GetHistoryDepth() << " name "
152  << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":" << module
153  << ":" << cell;
154 #endif
155  }
156 #ifdef EDM_ML_DEBUG
157  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
158  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
159  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
160  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
161  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
162  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
163  << touch->GetReplicaNumber(4) << " "
164  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
165  << mat->GetName() << ":" << mat->GetRadlen();
166 #endif
167  // The following statement should be examined later before elimination
168  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
169  return 0;
170 
171  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
172  if (!isItinFidVolume(hitPoint)) {
173 #ifdef EDM_ML_DEBUG
174  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCScintillatorDetId(id)
175  << " is rejected by fiducilal volume cut";
176 #endif
177  id = 0;
178  }
179  return id;
180 }
181 
182 void HGCScintSD::update(const BeginOfJob* job) {
183  if (hgcons_ != nullptr) {
186  levelT1_ = hgcons_->levelTop(0);
187  levelT2_ = hgcons_->levelTop(1);
188 #ifdef EDM_ML_DEBUG
189  edm::LogVerbatim("HGCSim") << "HGCScintSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
190  << " top Level " << levelT1_ << ":" << levelT2_;
191 #endif
192 
193  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_);
194  } else {
195  throw cms::Exception("Unknown", "HGCScintSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
196  }
197 }
198 
200 
201 bool HGCScintSD::filterHit(CaloG4Hit* aHit, double time) {
202  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
203 }
204 
205 uint32_t HGCScintSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
206  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
207  return id;
208 }
209 
210 bool HGCScintSD::isItinFidVolume(const G4ThreeVector& pos) {
211  if (fiducialCut_) {
212  return (hgcons_->distFromEdgeTrap(pos.x(), pos.y(), pos.z()) > distanceFromEdge_);
213  } else {
214  return true;
215  }
216 }
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
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:30
uint16_t *__restrict__ id
Definition: CaloSD.h:40
double distanceFromEdge_
Definition: HGCScintSD.h:46
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCScintSD.cc:123
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
const double MeV
double distFromEdgeTrap(double x, double y, double z) const
double weight_
Definition: HGCScintSD.h:50
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCScintSD.cc:210
double minSlope() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HGCalGeometryMode::GeometryMode geomMode() const
void initRun() override
Definition: HGCScintSD.cc:199
double tmaxHit
Definition: CaloSD.h:144
double getEnergyDeposit(const G4Step *) override
Definition: HGCScintSD.cc:90
double eminHit_
Definition: HGCScintSD.h:46
double slopeMin_
Definition: HGCScintSD.h:46
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCScintSD.h:42
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:657
std::string nameX_
Definition: HGCScintSD.h:44
void printDetectorLevels(const G4VTouchable *) const
Definition: CaloSD.cc:1072
DetId::Detector mydet_
Definition: HGCScintSD.h:43
int levelTop(int ind=0) const
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:849
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCScintSD.cc:182
tuple module
Definition: callgraph.py:69
double birk3_
Definition: HGCScintSD.h:50
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
bool storeAllG4Hits_
Definition: HGCScintSD.h:48
bool filterHit(CaloG4Hit *, double) override
Definition: HGCScintSD.cc:201