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