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 
23 #include "G4LogicalVolumeStore.hh"
24 #include "G4LogicalVolume.hh"
25 #include "G4Step.hh"
26 #include "G4Track.hh"
27 #include "G4ParticleTable.hh"
28 #include "G4VProcess.hh"
29 #include "G4Trap.hh"
30 
31 #include <iostream>
32 #include <fstream>
33 #include <iomanip>
34 
35 //#define EDM_ML_DEBUG
36 
38  const DDCompactView& cpv,
39  const SensitiveDetectorCatalog& clg,
40  edm::ParameterSet const& p,
41  const SimTrackManager* manager)
42  : CaloSD(name,
43  cpv,
44  clg,
45  p,
46  manager,
47  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
48  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
49  hgcons_(nullptr),
50  slopeMin_(0),
51  levelT1_(99),
52  levelT2_(99) {
53  numberingScheme_.reset(nullptr);
54 
55  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCScintSD");
56  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
57  fiducialCut_ = m_HGC.getParameter<bool>("FiducialCut");
58  distanceFromEdge_ = m_HGC.getParameter<double>("DistanceFromEdge");
59  useBirk_ = m_HGC.getParameter<bool>("UseBirkLaw");
60  birk1_ = m_HGC.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
61  birk2_ = m_HGC.getParameter<double>("BirkC2");
62  birk3_ = m_HGC.getParameter<double>("BirkC3");
63  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
64 
65  if (storeAllG4Hits_) {
66  setUseMap(false);
68  }
69 
70  //this is defined in the hgcsens.xml
71  G4String myName = name;
73  nameX_ = "HGCal";
74  if (myName.find("HitsHEback") != std::string::npos) {
76  nameX_ = "HGCalHEScintillatorSensitive";
77  }
78 
79 #ifdef EDM_ML_DEBUG
80  edm::LogVerbatim("HGCSim") << "**************************************************"
81  << "\n"
82  << "* *"
83  << "\n"
84  << "* Constructing a HGCScintSD with name " << name << "\n"
85  << "* *"
86  << "\n"
87  << "**************************************************";
88 #endif
89  edm::LogVerbatim("HGCSim") << "HGCScintSD:: Threshold for storing hits: " << eminHit_ << " for " << nameX_
90  << " detector " << mydet_;
91  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
92  edm::LogVerbatim("HGCSim") << "Fiducial volume cut with cut from eta/phi "
93  << "boundary " << fiducialCut_ << " at " << distanceFromEdge_;
94  edm::LogVerbatim("HGCSim") << "Use of Birks law is set to " << useBirk_
95  << " with three constants kB = " << birk1_ << ", C1 = " << birk2_ << ", C2 = " << birk3_;
96 }
97 
98 double HGCScintSD::getEnergyDeposit(const G4Step* aStep) {
99  double r = aStep->GetPreStepPoint()->GetPosition().perp();
100  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
101 #ifdef EDM_ML_DEBUG
102  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
103  G4String parName = aStep->GetTrack()->GetDefinition()->GetParticleName();
104  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
105  edm::LogVerbatim("HGCSim") << "HGCScintSD: Hit from standard path from " << lv->GetName() << " for Track "
106  << aStep->GetTrack()->GetTrackID() << " (" << parCode << ":" << parName << ") R = " << r
107  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
108 #endif
109  // Apply fiducial cut
110  if (r < z * slopeMin_) {
111 #ifdef EDM_ML_DEBUG
112  edm::LogVerbatim("HGCSim") << "HGCScintSD: Fiducial Volume cut";
113 #endif
114  return 0.0;
115  }
116 
117  double wt1 = getResponseWt(aStep->GetTrack());
118  double wt2 = aStep->GetTrack()->GetWeight();
119  double wt3 = (useBirk_ ? getAttenuation(aStep, birk1_, birk2_, birk3_) : 1.0);
120  double destep = weight_ * wt1 * wt3 * (aStep->GetTotalEnergyDeposit());
121  if (wt2 > 0)
122  destep *= wt2;
123 #ifdef EDM_ML_DEBUG
124  edm::LogVerbatim("HGCalSim") << "HGCScintSD: weights= " << weight_ << ":" << wt1 << ":" << wt2 << ":" << wt3
125  << " Total weight " << weight_ * wt1 * wt2 * wt3
126  << " deStep: " << aStep->GetTotalEnergyDeposit() << ":" << destep;
127 #endif
128  return destep;
129 }
130 
131 uint32_t HGCScintSD::setDetUnitId(const G4Step* aStep) {
132  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
133  const G4VTouchable* touch = preStepPoint->GetTouchable();
134 
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, module, cell;
141  if ((touch->GetHistoryDepth() == levelT1_) || (touch->GetHistoryDepth() == levelT2_)) {
142  layer = touch->GetReplicaNumber(0);
143  module = -1;
144  cell = -1;
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 
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 }
double birk2_
Definition: HGCScintSD.h:51
T getParameter(std::string const &) const
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCScintSD.h:46
int levelT2_
Definition: HGCScintSD.h:48
void setNumberCheckedHits(int val)
Definition: CaloSD.h:112
double birk1_
Definition: HGCScintSD.h:51
Definition: CaloSD.h:38
double distanceFromEdge_
Definition: HGCScintSD.h:47
#define nullptr
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCScintSD.cc:131
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
int levelT1_
Definition: HGCScintSD.h:48
void setUseMap(bool val)
Definition: CaloSD.h:101
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:49
bool useBirk_
Definition: HGCScintSD.h:50
const HGCalDDDConstants * hgcons_
Definition: HGCScintSD.h:42
HGCScintSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCScintSD.cc:37
const double MeV
double distFromEdgeTrap(double x, double y, double z) const
double weight_
Definition: HGCScintSD.h:51
bool isItinFidVolume(const G4ThreeVector &)
Definition: HGCScintSD.cc:218
double minSlope() const
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
HGCalGeometryMode::GeometryMode geomMode() const
void initRun() override
Definition: HGCScintSD.cc:207
double tmaxHit
Definition: CaloSD.h:133
double getEnergyDeposit(const G4Step *) override
Definition: HGCScintSD.cc:98
double eminHit_
Definition: HGCScintSD.h:47
double slopeMin_
Definition: HGCScintSD.h:47
std::unique_ptr< HGCalNumberingScheme > numberingScheme_
Definition: HGCScintSD.h:43
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:438
std::string nameX_
Definition: HGCScintSD.h:45
HLT enums.
T get() const
Definition: EventSetup.h:71
DetId::Detector mydet_
Definition: HGCScintSD.h:44
int levelTop(int ind=0) const
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:623
Definition: vlib.h:208
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCScintSD.cc:186
double birk3_
Definition: HGCScintSD.h:51
double getEnergyDeposit() const
Definition: CaloG4Hit.h:77
bool storeAllG4Hits_
Definition: HGCScintSD.h:49
bool filterHit(CaloG4Hit *, double) override
Definition: HGCScintSD.cc:209