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 
18 #include "G4LogicalVolumeStore.hh"
19 #include "G4LogicalVolume.hh"
20 #include "G4Step.hh"
21 #include "G4Track.hh"
22 #include "G4ParticleTable.hh"
23 #include "G4VProcess.hh"
24 #include "G4Trap.hh"
25 
26 #include <fstream>
27 #include <iomanip>
28 #include <iostream>
29 #include <memory>
30 
31 //#define EDM_ML_DEBUG
32 
34  const edm::EventSetup& es,
35  const SensitiveDetectorCatalog& clg,
36  edm::ParameterSet const& p,
37  const SimTrackManager* manager)
38  : CaloSD(name,
39  es,
40  clg,
41  p,
42  manager,
43  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
44  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
45  hgcons_(nullptr),
46  slopeMin_(0),
47  levelT1_(99),
48  levelT2_(99) {
49  numberingScheme_.reset(nullptr);
50 
51  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCScintSD");
52  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
53  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
54  distanceFromEdge_ = m_HGC.getParameter<double>("DistanceFromEdge");
55  useBirk_ = m_HGC.getParameter<bool>("UseBirkLaw");
56  birk1_ = m_HGC.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
57  birk2_ = m_HGC.getParameter<double>("BirkC2");
58  birk3_ = m_HGC.getParameter<double>("BirkC3");
59  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
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_;
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 
94 double HGCScintSD::getEnergyDeposit(const G4Step* aStep) {
95  double r = aStep->GetPreStepPoint()->GetPosition().perp();
96  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
97 #ifdef EDM_ML_DEBUG
98  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
99  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
100  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
101  edm::LogVerbatim("HGCSim") << "HGCScintSD: Hit from standard path from " << lv->GetName() << " for Track "
102  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
103  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
104 #endif
105  // Apply fiducial cut
106  if (r < z * slopeMin_) {
107 #ifdef EDM_ML_DEBUG
108  edm::LogVerbatim("HGCSim") << "HGCScintSD: Fiducial Volume cut";
109 #endif
110  return 0.0;
111  }
112 
113  double wt1 = getResponseWt(aStep->GetTrack());
114  double wt2 = aStep->GetTrack()->GetWeight();
115  double wt3 = (useBirk_ ? getAttenuation(aStep, birk1_, birk2_, birk3_) : 1.0);
116  double destep = weight_ * wt1 * wt3 * (aStep->GetTotalEnergyDeposit());
117  if (wt2 > 0)
118  destep *= wt2;
119 #ifdef EDM_ML_DEBUG
120  edm::LogVerbatim("HGCalSim") << "HGCScintSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << ":" << wt3
121  << " Total weight " << weight_ * wt1 * wt2 * wt3
122  << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
123 #endif
124  return destep;
125 }
126 
127 uint32_t HGCScintSD::setDetUnitId(const G4Step* aStep) {
128  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
129  const G4VTouchable* touch = preStepPoint->GetTouchable();
130 
131 #ifdef EDM_ML_DEBUG
132  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_;
133  printDetectorLevels(touch);
134 #endif
135  //determine the exact position in global coordinates in the mass geometry
136  G4ThreeVector hitPoint = preStepPoint->GetPosition();
137  float globalZ = touch->GetTranslation(0).z();
138  int iz(globalZ > 0 ? 1 : -1);
139 
140  int layer(0), module(-1), cell(-1);
142  layer = touch->GetReplicaNumber(1);
143  } else if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
144  layer = touch->GetReplicaNumber(0);
145 #ifdef EDM_ML_DEBUG
146  edm::LogVerbatim("HGCSim") << "DepthsTop: " << touch->GetHistoryDepth() << ":" << levelT1_ << ":" << levelT2_
147  << " name " << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":"
148  << module << ":" << cell;
149 #endif
150  } else {
151  layer = touch->GetReplicaNumber(3);
152  module = touch->GetReplicaNumber(2);
153  cell = touch->GetReplicaNumber(1);
154 #ifdef EDM_ML_DEBUG
155  edm::LogVerbatim("HGCSim") << "DepthsInside: " << touch->GetHistoryDepth() << " name "
156  << touch->GetVolume(0)->GetName() << " layer:module:cell " << layer << ":" << module
157  << ":" << cell;
158 #endif
159  }
160 #ifdef EDM_ML_DEBUG
161  G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
162  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
163  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
164  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
165  << touch->GetReplicaNumber(2) << " " << touch->GetVolume(3)->GetName() << ":"
166  << touch->GetReplicaNumber(3) << " " << touch->GetVolume(4)->GetName() << ":"
167  << touch->GetReplicaNumber(4) << " "
168  << " layer:module:cell " << layer << ":" << module << ":" << cell << " Material "
169  << mat->GetName() << ":" << mat->GetRadlen();
170 #endif
171  // The following statement should be examined later before elimination
172  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
173  return 0;
174 
175  uint32_t id = setDetUnitId(layer, module, cell, iz, hitPoint);
176  if (!isItinFidVolume(hitPoint)) {
177 #ifdef EDM_ML_DEBUG
178  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " " << HGCScintillatorDetId(id)
179  << " is rejected by fiducilal volume cut";
180 #endif
181  id = 0;
182  }
183  return id;
184 }
185 
186 void HGCScintSD::update(const BeginOfJob* job) {
187  const edm::EventSetup* es = (*job)();
189  es->get<IdealGeometryRecord>().get(nameX_, hdc);
190  if (hdc.isValid()) {
191  hgcons_ = hdc.product();
194  levelT1_ = hgcons_->levelTop(0);
195  levelT2_ = hgcons_->levelTop(1);
196 #ifdef EDM_ML_DEBUG
197  edm::LogVerbatim("HGCSim") << "HGCScintSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
198  << " top Level " << levelT1_ << ":" << levelT2_;
199 #endif
200 
201  numberingScheme_ = std::make_unique<HGCalNumberingScheme>(*hgcons_, mydet_, nameX_);
202  } else {
203  throw cms::Exception("Unknown", "HGCScintSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
204  }
205 }
206 
208 
209 bool HGCScintSD::filterHit(CaloG4Hit* aHit, double time) {
210  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
211 }
212 
213 uint32_t HGCScintSD::setDetUnitId(int layer, int module, int cell, int iz, G4ThreeVector& pos) {
214  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(layer, module, cell, iz, pos, weight_) : 0;
215  return id;
216 }
217 
218 bool HGCScintSD::isItinFidVolume(const G4ThreeVector& pos) {
219  if (fiducialCut_) {
220  return (hgcons_->distFromEdgeTrap(pos.x(), pos.y(), pos.z()) > distanceFromEdge_);
221  } else {
222  return true;
223  }
224 }
HGCScintSD::eminHit_
double eminHit_
Definition: HGCScintSD.h:46
SimTrackManager
Definition: SimTrackManager.h:35
electrons_cff.bool
bool
Definition: electrons_cff.py:366
CaloSD::tmaxHit
double tmaxHit
Definition: CaloSD.h:141
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
HGCScintillatorDetId.h
ESHandle.h
HGCScintSD::weight_
double weight_
Definition: HGCScintSD.h:50
HGCScintSD::isItinFidVolume
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCScintSD.cc:218
HGCalDDDConstants::geomMode
HGCalGeometryMode::GeometryMode geomMode() const
Definition: HGCalDDDConstants.h:54
edm
HLT enums.
Definition: AlignableModifier.h:19
HGCalGeometryMode.h
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
pos
Definition: PixelAliasList.h:18
HGCalDDDConstants::distFromEdgeTrap
double distFromEdgeTrap(double x, double y, double z) const
Definition: HGCalDDDConstants.cc:325
CaloSD::getAttenuation
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:685
protons_cff.time
time
Definition: protons_cff.py:39
MeV
const double MeV
HGCScintSD::slopeMin_
double slopeMin_
Definition: HGCScintSD.h:46
HGCScintSD::getEnergyDeposit
double getEnergyDeposit(const G4Step *) override
Definition: HGCScintSD.cc:94
FastMath.h
CaloSD::getResponseWt
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:876
edm::EventSetup::get
T get() const
Definition: EventSetup.h:87
HGCScintSD::numberingScheme_
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCScintSD.h:42
CaloG4Hit::getEnergyDeposit
double getEnergyDeposit() const
Definition: CaloG4Hit.h:78
DDAxes::z
HGCScintSD::mydet_
DetId::Detector mydet_
Definition: HGCScintSD.h:43
edm::ESHandle
Definition: DTSurvey.h:22
HGCScintSD::filterHit
bool filterHit(CaloG4Hit *, double) override
Definition: HGCScintSD.cc:209
HGCScintSD::birk2_
double birk2_
Definition: HGCScintSD.h:50
HGCScintSD::nameX_
std::string nameX_
Definition: HGCScintSD.h:44
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
BeginOfJob
Definition: BeginOfJob.h:8
HGCScintSD::geom_mode_
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCScintSD.h:45
HGCScintSD::update
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCScintSD.cc:186
phase1PixelTopology::layer
constexpr std::array< uint8_t, layerIndexSize > layer
Definition: phase1PixelTopology.h:99
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HGCScintSD::birk3_
double birk3_
Definition: HGCScintSD.h:50
HGCScintSD::birk1_
double birk1_
Definition: HGCScintSD.h:50
HGCScintSD::storeAllG4Hits_
bool storeAllG4Hits_
Definition: HGCScintSD.h:48
edm::ParameterSet
Definition: ParameterSet.h:47
HGCalGeometryMode::TrapezoidModule
Definition: HGCalGeometryMode.h:35
ParameterSet
Definition: Functions.h:16
HGCScintSD::distanceFromEdge_
double distanceFromEdge_
Definition: HGCScintSD.h:46
CaloSD::printDetectorLevels
void printDetectorLevels(const G4VTouchable *) const
Definition: CaloSD.cc:1104
CaloSD::setNumberCheckedHits
void setNumberCheckedHits(int val)
Definition: CaloSD.h:119
IdealGeometryRecord.h
HGCScintSD::setDetUnitId
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCScintSD.cc:127
edm::EventSetup
Definition: EventSetup.h:58
TrackInformation.h
CaloG4Hit
Definition: CaloG4Hit.h:32
get
#define get
HGCScintSD::levelT2_
int levelT2_
Definition: HGCScintSD.h:47
HGCScintSD::HGCScintSD
HGCScintSD(const std::string &, const edm::EventSetup &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCScintSD.cc:33
alignCSCRings.r
r
Definition: alignCSCRings.py:93
HGCalDDDConstants::minSlope
double minSlope() const
Definition: HGCalDDDConstants.h:96
HGCScintSD::useBirk_
bool useBirk_
Definition: HGCScintSD.h:49
HGCScintSD::levelT1_
int levelT1_
Definition: HGCScintSD.h:47
HGCScintillatorDetId
Definition: HGCScintillatorDetId.h:23
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:29
DetId::HGCalHSc
Definition: DetId.h:34
HGCScintSD::hgcons_
const HGCalDDDConstants * hgcons_
Definition: HGCScintSD.h:41
Exception
Definition: hltDiff.cc:245
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
EventSetup.h
HGCScintSD.h
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
Exception.h
HGCalDDDConstants::levelTop
int levelTop(int ind=0) const
Definition: HGCalDDDConstants.h:88
CaloSD::setUseMap
void setUseMap(bool val)
Definition: CaloSD.h:108
funct::abs
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HGCalDDDConstants.h
HGCScintSD::initRun
void initRun() override
Definition: HGCScintSD.cc:207
DetId::Forward
Definition: DetId.h:30
TauDecayModes.dec
dec
Definition: TauDecayModes.py:143
HGCScintSD::fiducialCut_
bool fiducialCut_
Definition: HGCScintSD.h:48
IdealGeometryRecord
Definition: IdealGeometryRecord.h:25
CaloSD
Definition: CaloSD.h:38
g
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
HGCalTestNumbering.h