CMS 3D CMS Logo

HGCSD.cc
Go to the documentation of this file.
1 // File: HGCSD.cc
3 // Description: Sensitive Detector class for Combined Forward Calorimeter
5 
7 
16 #include "G4LogicalVolumeStore.hh"
17 #include "G4LogicalVolume.hh"
18 #include "G4Step.hh"
19 #include "G4Track.hh"
20 #include "G4ParticleTable.hh"
21 #include "G4VProcess.hh"
22 #include "G4Trap.hh"
23 
24 #include <fstream>
25 #include <iomanip>
26 #include <iostream>
27 #include <memory>
28 
29 //#define EDM_ML_DEBUG
30 //#define plotDebug
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  levelT_(99),
46  tree_(nullptr) {
47  numberingScheme_.reset(nullptr);
48  mouseBite_.reset(nullptr);
49 
50  edm::ParameterSet m_HGC = p.getParameter<edm::ParameterSet>("HGCSD");
51  eminHit_ = m_HGC.getParameter<double>("EminHit") * CLHEP::MeV;
52  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
53  rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
54  waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
55  angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
56  double waferSize = m_HGC.getUntrackedParameter<double>("WaferSize") * CLHEP::mm;
57  double mouseBite = m_HGC.getUntrackedParameter<double>("MouseBite") * CLHEP::mm;
58  mouseBiteCut_ = waferSize * tan(30.0 * CLHEP::deg) - mouseBite;
59  dd4hep_ = p.getParameter<bool>("g4GeometryDD4hepSource");
60 
61  if (storeAllG4Hits_) {
62  setUseMap(false);
64  }
65  //this is defined in the hgcsens.xml
66  G4String myName = name;
68  nameX_ = "HGCal";
69  if (myName.find("HitsEE") != std::string::npos) {
71  nameX_ = "HGCalEESensitive";
72  } else if (myName.find("HitsHEfront") != std::string::npos) {
74  nameX_ = "HGCalHESiliconSensitive";
75  } else if (myName.find("HitsHEback") != std::string::npos) {
77  nameX_ = "HGCalHEScintillatorSensitive";
78  }
79 
80 #ifdef EDM_ML_DEBUG
81  edm::LogVerbatim("HGCSim") << "**************************************************"
82  << "\n"
83  << "* *"
84  << "\n"
85  << "* Constructing a HGCSD with name " << name << "\n"
86  << "* *"
87  << "\n"
88  << "**************************************************";
89 #endif
90  edm::LogVerbatim("HGCSim") << "HGCSD:: Threshold for storing hits: " << eminHit << " for " << nameX_ << " subdet "
91  << myFwdSubdet_;
92  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits " << storeAllG4Hits_;
93  edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_ << " Size of wafer " << waferSize
94  << " Mouse Bite " << mouseBite << ":" << mouseBiteCut_ << " along " << angles_.size()
95  << " axes";
96 }
97 
98 double HGCSD::getEnergyDeposit(const G4Step* aStep) {
99  double r = aStep->GetPreStepPoint()->GetPosition().perp();
100  double z = std::abs(aStep->GetPreStepPoint()->GetPosition().z());
101 
102 #ifdef EDM_ML_DEBUG
103  G4int parCode = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
104  G4LogicalVolume* lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
105  edm::LogVerbatim("HGCSim") << "HGCSD: Hit from standard path from " << lv->GetName() << " for Track "
106  << aStep->GetTrack()->GetTrackID() << " ("
107  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ":" << parCode << ") R = " << r
108  << " Z = " << z << " slope = " << r / z << ":" << slopeMin_;
109 #endif
110 
111  // Apply fiductial volume
112  if (r < z * slopeMin_) {
113  return 0.0;
114  }
115 
116  double wt1 = getResponseWt(aStep->GetTrack());
117  double wt2 = aStep->GetTrack()->GetWeight();
118  double destep = wt1 * aStep->GetTotalEnergyDeposit();
119  if (wt2 > 0)
120  destep *= wt2;
121 
122 #ifdef plotDebug
123  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
124  G4double tmptrackE = aStep->GetTrack()->GetKineticEnergy();
125  G4int parCodex = aStep->GetTrack()->GetDefinition()->GetPDGEncoding();
126  G4double angle = (aStep->GetTrack()->GetMomentumDirection().theta()) / CLHEP::deg;
127  G4int layer = ((touch->GetHistoryDepth() == levelT_) ? touch->GetReplicaNumber(0) : touch->GetReplicaNumber(2));
128  G4int ilayer = (layer - 1) / 3;
129  if (aStep->GetTotalEnergyDeposit() > 0) {
130  t_Layer_.emplace_back(ilayer);
131  t_Parcode_.emplace_back(parCodex);
132  t_dEStep1_.emplace_back(aStep->GetTotalEnergyDeposit());
133  t_dEStep2_.emplace_back(destep);
134  t_TrackE_.emplace_back(tmptrackE);
135  t_Angle_.emplace_back(angle);
136  }
137 #endif
138 
139  return destep;
140 }
141 
142 uint32_t HGCSD::setDetUnitId(const G4Step* aStep) {
143  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
144  const G4VTouchable* touch = preStepPoint->GetTouchable();
145 
146  //determine the exact position in global coordinates in the mass geometry
147  G4ThreeVector hitPoint = preStepPoint->GetPosition();
148  float globalZ = touch->GetTranslation(0).z();
149  int iz(globalZ > 0 ? 1 : -1);
150 
151  //convert to local coordinates (=local to the current volume):
152  G4ThreeVector localpos = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
153 
154  //get the det unit id with
156 
157  int layer(-1), moduleLev(-1), module(-1), cell(-1);
158  if (touch->GetHistoryDepth() == levelT_) {
159  layer = touch->GetReplicaNumber(0);
160 #ifdef EDM_ML_DEBUG
161  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
162  << " layer:module:cell " << layer << ":" << moduleLev << ":" << module << ":" << cell;
163 #endif
164  } else {
165  layer = touch->GetReplicaNumber(2);
166  module = touch->GetReplicaNumber(1);
167  cell = touch->GetReplicaNumber(0);
168  moduleLev = 1;
169  }
170 #ifdef EDM_ML_DEBUG
171  const G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
172  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth() << " name " << touch->GetVolume(0)->GetName()
173  << ":" << touch->GetReplicaNumber(0) << " " << touch->GetVolume(1)->GetName() << ":"
174  << touch->GetReplicaNumber(1) << " " << touch->GetVolume(2)->GetName() << ":"
175  << touch->GetReplicaNumber(2) << " layer:module:cell " << layer << ":" << moduleLev
176  << ":" << module << ":" << cell << " Material " << mat->GetName() << ":"
177  << mat->GetRadlen();
178  for (int k = 0; k < touch->GetHistoryDepth(); ++k)
179  edm::LogVerbatim("HGCSim") << "Level [" << k << "] " << touch->GetVolume(k)->GetName() << ":"
180  << touch->GetReplicaNumber(k);
181 #endif
182  // The following statement should be examined later before elimination
183  // VI: this is likely a check if media is vacuum - not needed
184  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.)
185  return 0;
186 
187  uint32_t id = setDetUnitId(subdet, layer, module, cell, iz, localpos);
188  if (rejectMB_ && id != 0) {
189  int det, z, lay, wafer, type, ic;
190  HGCalTestNumbering::unpackHexagonIndex(id, det, z, lay, wafer, type, ic);
191 #ifdef EDM_ML_DEBUG
192  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec << " Input " << subdet << ":" << layer << ":"
193  << module << ":" << cell << ":" << iz << localpos.x() << ":" << localpos.y()
194  << " Decode " << det << ":" << z << ":" << lay << ":" << wafer << ":" << type << ":"
195  << ic;
196 #endif
197  G4ThreeVector local =
198  ((moduleLev >= 0) ? (touch->GetHistory()->GetTransform(moduleLev).TransformPoint(hitPoint)) : G4ThreeVector());
199  if (mouseBite_->exclude(local, z, layer, wafer, 0))
200  id = 0;
201  }
202  return id;
203 }
204 
205 void HGCSD::update(const BeginOfJob* job) {
206  if (hgcons_ != nullptr) {
209  levelT_ = hgcons_->levelTop();
210  if (dd4hep_)
211  ++levelT_;
212  numberingScheme_ = std::make_unique<HGCNumberingScheme>(*hgcons_, nameX_);
213  if (rejectMB_)
214  mouseBite_ = std::make_unique<HGCMouseBite>(*hgcons_, angles_, mouseBiteCut_, waferRot_);
215  } else {
216  edm::LogError("HGCSim") << "HGCSD : Cannot find HGCalDDDConstants for " << nameX_;
217  throw cms::Exception("Unknown", "HGCSD") << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
218  }
219 #ifdef EDM_ML_DEBUG
220  edm::LogVerbatim("HGCSim") << "HGCSD::Initialized with mode " << geom_mode_ << " Slope cut " << slopeMin_
221  << " top Level " << levelT_ << " dd4hep " << dd4hep_;
222 #endif
223 }
224 
226 #ifdef plotDebug
228  if (tfile.isAvailable()) {
229  tree_ = tfile->make<TTree>("TreeHGCSD", "TreeHGCSD");
230  tree_->Branch("EventID", &t_EventID_);
231  tree_->Branch("Layer", &t_Layer_);
232  tree_->Branch("ParticleCode", &t_Parcode_);
233  tree_->Branch("dEStepOriginal", &t_dEStep1_);
234  tree_->Branch("dEStepWeighted", &t_dEStep2_);
235  tree_->Branch("TrackEnergy", &t_TrackE_);
236  tree_->Branch("ThetaAngle", &t_Angle_);
237  }
238 #endif
239 }
240 
241 void HGCSD::initEvent(const BeginOfEvent* g4Event) {
242  const G4Event* evt = (*g4Event)();
243  t_EventID_ = evt->GetEventID();
244 #ifdef plotDebug
245  t_Layer_.clear();
246  t_Parcode_.clear();
247  t_dEStep1_.clear();
248  t_dEStep2_.clear();
249  t_TrackE_.clear();
250  t_Angle_.clear();
251 #endif
252 }
253 
255 #ifdef plotDebug
256  if (tree_)
257  tree_->Fill();
258 #endif
259 }
260 
261 bool HGCSD::filterHit(CaloG4Hit* aHit, double time) {
262  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
263 }
264 
265 uint32_t HGCSD::setDetUnitId(ForwardSubdetector& subdet, int layer, int module, int cell, int iz, G4ThreeVector& pos) {
266  uint32_t id = numberingScheme_ ? numberingScheme_->getUnitID(subdet, layer, module, cell, iz, pos) : 0;
267  return id;
268 }
std::vector< double > t_TrackE_
Definition: HGCSD.h:63
bool dd4hep_
Definition: HGCSD.h:57
Log< level::Info, true > LogVerbatim
double getEnergyDeposit(const G4Step *) override
Definition: HGCSD.cc:98
std::vector< double > t_dEStep1_
Definition: HGCSD.h:63
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
void setNumberCheckedHits(int val)
Definition: CaloSD.h:122
std::vector< double > t_dEStep2_
Definition: HGCSD.h:63
int levelT_
Definition: HGCSD.h:54
double mouseBiteCut_
Definition: HGCSD.h:56
Definition: CaloSD.h:40
std::vector< int > t_Layer_
Definition: HGCSD.h:62
std::string nameX_
Definition: HGCSD.h:47
std::vector< double > angles_
Definition: HGCSD.h:58
Log< level::Error, false > LogError
bool storeAllG4Hits_
Definition: HGCSD.h:55
HGCalGeometryMode::GeometryMode geomMode() const
ForwardSubdetector
bool waferRot_
Definition: HGCSD.h:55
void setUseMap(bool val)
Definition: CaloSD.h:111
double eminHit
Definition: CaloSD.h:144
T getUntrackedParameter(std::string const &, T const &) const
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCSD.h:48
double slopeMin_
Definition: HGCSD.h:53
std::unique_ptr< HGCMouseBite > mouseBite_
Definition: HGCSD.h:50
Definition: tfile.py:1
const HGCalDDDConstants * hgcons_
Definition: HGCSD.h:46
void initRun() override
Definition: HGCSD.cc:225
Tan< T >::type tan(const T &t)
Definition: Tan.h:22
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< int > t_Parcode_
Definition: HGCSD.h:62
std::unique_ptr< HGCNumberingScheme > numberingScheme_
Definition: HGCSD.h:49
void initEvent(const BeginOfEvent *) override
Definition: HGCSD.cc:241
double tmaxHit
Definition: CaloSD.h:144
double eminHit_
Definition: HGCSD.h:51
std::vector< double > t_Angle_
Definition: HGCSD.h:64
bool rejectMB_
Definition: HGCSD.h:55
double minSlope() const
TTree * tree_
Definition: HGCSD.h:60
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
HGCSD(const std::string &, const HGCalDDDConstants *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCSD.cc:32
HLT enums.
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCSD.cc:142
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCSD.cc:205
bool filterHit(CaloG4Hit *, double) override
Definition: HGCSD.cc:261
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:860
int levelTop(int ind=0) const
ForwardSubdetector myFwdSubdet_
Definition: HGCSD.h:52
void endEvent() override
Definition: HGCSD.cc:254
static void unpackHexagonIndex(const uint32_t &idx, int &subdet, int &z, int &lay, int &wafer, int &celltyp, int &cell)
uint32_t t_EventID_
Definition: HGCSD.h:61
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11