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 
22 #include "G4LogicalVolumeStore.hh"
23 #include "G4LogicalVolume.hh"
24 #include "G4Step.hh"
25 #include "G4Track.hh"
26 #include "G4ParticleTable.hh"
27 #include "G4VProcess.hh"
28 #include "G4Trap.hh"
29 
30 #include <iostream>
31 #include <fstream>
32 #include <iomanip>
33 
34 //#define EDM_ML_DEBUG
35 
37  const SensitiveDetectorCatalog & clg,
38  edm::ParameterSet const & p, const SimTrackManager* manager) :
39  CaloSD(name, cpv, clg, p, manager,
40  (float)(p.getParameter<edm::ParameterSet>("HGCSD").getParameter<double>("TimeSliceUnit")),
41  p.getParameter<edm::ParameterSet>("HGCSD").getParameter<bool>("IgnoreTrackID")),
42  slopeMin_(0), levelT_(99) {
43 
44  numberingScheme_.reset(nullptr); mouseBite_.reset(nullptr);
45 
47  eminHit_ = m_HGC.getParameter<double>("EminHit")*CLHEP::MeV;
48  storeAllG4Hits_ = m_HGC.getParameter<bool>("StoreAllG4Hits");
49  rejectMB_ = m_HGC.getParameter<bool>("RejectMouseBite");
50  waferRot_ = m_HGC.getParameter<bool>("RotatedWafer");
51  angles_ = m_HGC.getUntrackedParameter<std::vector<double>>("WaferAngles");
52  double waferSize = m_HGC.getUntrackedParameter<double>("WaferSize")*CLHEP::mm;
53  double mouseBite = m_HGC.getUntrackedParameter<double>("MouseBite")*CLHEP::mm;
54  mouseBiteCut_ = waferSize*tan(30.0*CLHEP::deg) - mouseBite;
55 
56  if(storeAllG4Hits_) {
57  setUseMap(false);
59  }
60  //this is defined in the hgcsens.xml
61  G4String myName = name;
63  nameX_ = "HGCal";
64  if (myName.find("HitsEE")!=std::string::npos) {
66  nameX_ = "HGCalEESensitive";
67  } else if (myName.find("HitsHEfront")!=std::string::npos) {
69  nameX_ = "HGCalHESiliconSensitive";
70  } else 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 HGCSD with name " << name << "\n"
81  << "* *"
82  << "\n"
83  << "**************************************************";
84 #endif
85  edm::LogVerbatim("HGCSim") << "HGCSD:: Threshold for storing hits: "
86  << eminHit << " for " << nameX_ << " subdet "
87  << myFwdSubdet_;
88  edm::LogVerbatim("HGCSim") << "Flag for storing individual Geant4 Hits "
89  << storeAllG4Hits_;
90  edm::LogVerbatim("HGCSim") << "Reject MosueBite Flag: " << rejectMB_
91  << " Size of wafer " << waferSize
92  << " Mouse Bite " << mouseBite << ":"
93  << mouseBiteCut_ << " along " << angles_.size()
94  << " axes";
95 }
96 
97 double HGCSD::getEnergyDeposit(const G4Step* aStep) {
98 
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 =
105  aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
106  edm::LogVerbatim("HGCSim") << "HGCSD: Hit from standard path from "
107  << lv->GetName() << " for Track "
108  << aStep->GetTrack()->GetTrackID() << " ("
109  << aStep->GetTrack()->GetDefinition()->GetParticleName()
110  << ") R = " << r << " Z = "
111  << z << " slope = " << r/z << ":" << slopeMin_;
112 #endif
113 
114  // Apply fiductial volume
115  if (r < z*slopeMin_) { return 0.0; }
116 
117  double wt1 = getResponseWt(aStep->GetTrack());
118  double wt2 = aStep->GetTrack()->GetWeight();
119  double destep = wt1*aStep->GetTotalEnergyDeposit();
120  if (wt2 > 0) destep *= wt2;
121  return destep;
122 }
123 
124 uint32_t HGCSD::setDetUnitId(const G4Step * aStep) {
125 
126  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
127  const G4VTouchable* touch = preStepPoint->GetTouchable();
128 
129  //determine the exact position in global coordinates in the mass geometry
130  G4ThreeVector hitPoint = preStepPoint->GetPosition();
131  float globalZ=touch->GetTranslation(0).z();
132  int iz( globalZ>0 ? 1 : -1);
133 
134  //convert to local coordinates (=local to the current volume):
135  G4ThreeVector localpos = touch->GetHistory()->GetTopTransform().TransformPoint(hitPoint);
136 
137  //get the det unit id with
139 
140  int layer, module, cell;
141  if (touch->GetHistoryDepth() == levelT_) {
142  layer = touch->GetReplicaNumber(0);
143  module = -1;
144  cell = -1;
145 #ifdef EDM_ML_DEBUG
146  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth()
147  << " name " << touch->GetVolume(0)->GetName()
148  << " layer:module:cell " << layer << ":"
149  << module << ":" << cell;
150 #endif
151  } else {
152  layer = touch->GetReplicaNumber(2);
153  module = touch->GetReplicaNumber(1);
154  cell = touch->GetReplicaNumber(0);
155  }
156 #ifdef EDM_ML_DEBUG
157  edm::LogVerbatim("HGCSim") << "Depths: " << touch->GetHistoryDepth()
158  << " name " << touch->GetVolume(0)->GetName()
159  << ":" << touch->GetReplicaNumber(0) << " "
160  << touch->GetVolume(1)->GetName()
161  << ":" << touch->GetReplicaNumber(1) << " "
162  << touch->GetVolume(2)->GetName()
163  << ":" << touch->GetReplicaNumber(2) << " "
164  << " layer:module:cell " << layer << ":"
165  << module << ":" << cell <<" Material "
166  << mat->GetName() << ":"
167  << aStep->GetPreStepPoint()->GetMaterial()->GetRadlen();
168 #endif
169  // The following statement should be examined later before elimination
170  // VI: this is likely a check if media is vacuum - not needed
171  if (aStep->GetPreStepPoint()->GetMaterial()->GetRadlen() > 100000.) return 0;
172 
173  uint32_t id = setDetUnitId (subdet, layer, module, cell, iz, localpos);
174  if (rejectMB_ && id != 0) {
175  int det, z, lay, wafer, type, ic;
176  HGCalTestNumbering::unpackHexagonIndex(id, det, z, lay, wafer, type, ic);
177 #ifdef EDM_ML_DEBUG
178  edm::LogVerbatim("HGCSim") << "ID " << std::hex << id << std::dec
179  << " Decode " << det << ":" << z << ":" << lay
180  << ":" << wafer << ":" << type << ":" << ic;
181 #endif
182  if (mouseBite_->exclude(hitPoint, z, wafer, 0)) id = 0;
183  }
184  return id;
185 }
186 
187 void HGCSD::update(const BeginOfJob * job) {
188 
189  const edm::EventSetup* es = (*job)();
191  es->get<IdealGeometryRecord>().get(nameX_,hdc);
192  if (hdc.isValid()) {
193  const HGCalDDDConstants* hgcons = hdc.product();
194  geom_mode_ = hgcons->geomMode();
195  slopeMin_ = hgcons->minSlope();
196  levelT_ = hgcons->levelTop();
197  numberingScheme_.reset(new HGCNumberingScheme(*hgcons,nameX_));
198  if (rejectMB_) mouseBite_.reset(new HGCMouseBite(*hgcons,angles_,
200  } else {
201  edm::LogError("HGCSim") << "HCalSD : Cannot find HGCalDDDConstants for "
202  << nameX_;
203  throw cms::Exception("Unknown", "HGCSD")
204  << "Cannot find HGCalDDDConstants for " << nameX_ << "\n";
205  }
206 #ifdef EDM_ML_DEBUG
207  edm::LogVerbatim("HGCSim") << "HGCSD::Initialized with mode " << geom_mode_
208  << " Slope cut " << slopeMin_ << " top Level "
209  << levelT_;
210 #endif
211 }
212 
214 }
215 
216 bool HGCSD::filterHit(CaloG4Hit* aHit, double time) {
217  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit_));
218 }
219 
220 uint32_t HGCSD::setDetUnitId (ForwardSubdetector &subdet, int layer, int module,
221  int cell, int iz, G4ThreeVector &pos) {
222  uint32_t id = numberingScheme_ ?
223  numberingScheme_->getUnitID(subdet, layer, module, cell, iz, pos) : 0;
224  return id;
225 }
226 
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:97
void setNumberCheckedHits(int val)
Definition: CaloSD.h:109
int levelT_
Definition: HGCSD.h:54
double mouseBiteCut_
Definition: HGCSD.h:56
Definition: CaloSD.h:37
std::string nameX_
Definition: HGCSD.h:47
std::vector< double > angles_
Definition: HGCSD.h:57
bool storeAllG4Hits_
Definition: HGCSD.h:55
ForwardSubdetector
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
bool waferRot_
Definition: HGCSD.h:55
void setUseMap(bool val)
Definition: CaloSD.h:98
double eminHit
Definition: CaloSD.h:132
const double MeV
HGCalGeometryMode::GeometryMode geom_mode_
Definition: HGCSD.h:48
double slopeMin_
Definition: HGCSD.h:53
std::unique_ptr< HGCMouseBite > mouseBite_
Definition: HGCSD.h:50
void initRun() override
Definition: HGCSD.cc:213
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::unique_ptr< HGCNumberingScheme > numberingScheme_
Definition: HGCSD.h:49
double tmaxHit
Definition: CaloSD.h:132
double eminHit_
Definition: HGCSD.h:51
bool rejectMB_
Definition: HGCSD.h:55
HGCSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: HGCSD.cc:36
HLT enums.
uint32_t setDetUnitId(const G4Step *step) override
Definition: HGCSD.cc:124
T get() const
Definition: EventSetup.h:63
void update(const BeginOfJob *) override
This routine will be called when the appropriate signal arrives.
Definition: HGCSD.cc:187
bool filterHit(CaloG4Hit *, double) override
Definition: HGCSD.cc:216
int levelTop(int ind=0) const
double getResponseWt(const G4Track *)
Definition: CaloSD.cc:646
ForwardSubdetector myFwdSubdet_
Definition: HGCSD.h:52
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:81