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