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