CMS 3D CMS Logo

AHCalSD.cc
Go to the documentation of this file.
5 
8 
9 #include "G4LogicalVolume.hh"
10 #include "G4LogicalVolumeStore.hh"
11 #include "G4ParticleTable.hh"
12 #include "G4Step.hh"
13 #include "G4Track.hh"
14 #include "G4VProcess.hh"
15 
16 #include "G4PhysicalConstants.hh"
17 #include "G4SystemOfUnits.hh"
18 
19 #include <iomanip>
20 #include <map>
21 #include <string>
22 
23 //#define EDM_ML_DEBUG
24 
25 class AHCalSD : public CaloSD {
26 public:
27  AHCalSD(const std::string&,
28  const edm::EventSetup&,
30  edm::ParameterSet const&,
31  const SimTrackManager*);
32  ~AHCalSD() override = default;
33  uint32_t setDetUnitId(const G4Step* step) override;
34  bool unpackIndex(const uint32_t& idx, int& row, int& col, int& depth);
35 
36 protected:
37  double getEnergyDeposit(const G4Step*) override;
38  bool filterHit(CaloG4Hit*, double) override;
39 
40 private:
41  bool useBirk;
42  double birk1, birk2, birk3, betaThr;
43  double eminHit;
44 };
45 
47  const edm::EventSetup& es,
48  const SensitiveDetectorCatalog& clg,
49  edm::ParameterSet const& p,
50  const SimTrackManager* manager)
51  : CaloSD(name,
52  es,
53  clg,
54  p,
55  manager,
56  (float)(p.getParameter<edm::ParameterSet>("AHCalSD").getParameter<double>("TimeSliceUnit")),
57  p.getParameter<edm::ParameterSet>("AHCalSD").getParameter<bool>("IgnoreTrackID")) {
58  edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("AHCalSD");
59  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
60  birk1 = m_HC.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
61  birk2 = m_HC.getParameter<double>("BirkC2");
62  birk3 = m_HC.getParameter<double>("BirkC3");
63  eminHit = m_HC.getParameter<double>("EminHit") * CLHEP::MeV;
64 
65  edm::LogVerbatim("HcalSim") << "AHCalSD:: Use of Birks law is set to " << useBirk
66  << " with three constants kB = " << birk1 << ", C1 = " << birk2 << ", C2 = " << birk3
67  << "\nAHCalSD:: Threshold for storing"
68  << " hits: " << eminHit;
69 }
70 
71 double AHCalSD::getEnergyDeposit(const G4Step* aStep) {
72  double destep = aStep->GetTotalEnergyDeposit();
73  double wt2 = aStep->GetTrack()->GetWeight();
74  double weight = (wt2 > 0.0) ? wt2 : 1.0;
75 #ifdef EDM_ML_DEBUG
76  double weight0 = weight;
77 #endif
78  if (useBirk)
79  weight *= getAttenuation(aStep, birk1, birk2, birk3);
80 #ifdef EDM_ML_DEBUG
81  edm::LogVerbatim("HcalSim") << "AHCalSD: weight " << weight0 << " " << weight;
82 #endif
83  double edep = weight * destep;
84  return edep;
85 }
86 
87 uint32_t AHCalSD::setDetUnitId(const G4Step* aStep) {
88  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
89  const G4VTouchable* touch = preStepPoint->GetTouchable();
90 
91  int depth = (touch->GetReplicaNumber(1));
92  int incol = ((touch->GetReplicaNumber(0)) % AHCalParameters::kColumn_);
93  int inrow = ((touch->GetReplicaNumber(0)) / AHCalParameters::kColumn_) % AHCalParameters::kRow_;
94  int jncol = ((touch->GetReplicaNumber(0)) / AHCalParameters::kRowColumn_) % AHCalParameters::kSign_;
95  int jnrow = ((touch->GetReplicaNumber(0)) / AHCalParameters::kSignRowColumn_) % AHCalParameters::kSign_;
96  int col = (jncol == 0) ? incol : -incol;
97  int row = (jnrow == 0) ? inrow : -inrow;
98  uint32_t index = AHCalDetId(row, col, depth).rawId();
99 #ifdef EDM_ML_DEBUG
100  edm::LogVerbatim("HcalSim") << "AHCalSD: det = " << HcalOther << " depth = " << depth << " row = " << row
101  << " column = " << col << " packed index = 0x" << std::hex << index << std::dec;
102  bool flag = unpackIndex(index, row, col, depth);
103  edm::LogVerbatim("HcalSim") << "Results from unpacker for 0x" << std::hex << index << std::dec << " Flag " << flag
104  << " Row " << row << " Col " << col << " Depth " << depth;
105 #endif
106  return index;
107 }
108 
109 bool AHCalSD::unpackIndex(const uint32_t& idx, int& row, int& col, int& depth) {
110  DetId gen(idx);
111  HcalSubdetector subdet = (HcalSubdetector(gen.subdetId()));
112  bool rcode = (gen.det() == DetId::Hcal && subdet != HcalOther);
113  row = col = depth = 0;
114  if (rcode) {
115  row = AHCalDetId(idx).irow();
116  col = AHCalDetId(idx).icol();
117  depth = AHCalDetId(idx).depth();
118  }
119 #ifdef EDM_ML_DEBUG
120  edm::LogVerbatim("HcalSim") << "AHCalSD: packed index = 0x" << std::hex << idx << std::dec << " Row " << row
121  << " Column " << col << " Depth " << depth << " OK " << rcode;
122 #endif
123  return rcode;
124 }
125 
126 bool AHCalSD::filterHit(CaloG4Hit* aHit, double time) {
127  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit));
128 }
129 
133 
HcalOther
Definition: HcalAssistant.h:38
SimTrackManager
Definition: SimTrackManager.h:35
electrons_cff.bool
bool
Definition: electrons_cff.py:393
CaloSD::tmaxHit
double tmaxHit
Definition: CaloSD.h:137
dqmMemoryStats.float
float
Definition: dqmMemoryStats.py:127
AHCalSD::birk2
double birk2
Definition: AHCalSD.cc:42
AHCalSD::~AHCalSD
~AHCalSD() override=default
AHCalSD::unpackIndex
bool unpackIndex(const uint32_t &idx, int &row, int &col, int &depth)
Definition: AHCalSD.cc:109
AHCalDetId::irow
int irow() const
get the row number
Definition: AHCalDetId.cc:29
step
step
Definition: StallMonitor.cc:94
edm
HLT enums.
Definition: AlignableModifier.h:19
mps_merge.weight
weight
Definition: mps_merge.py:88
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
cuy.col
col
Definition: cuy.py:1010
DetId::Hcal
Definition: DetId.h:28
AHCalSD::setDetUnitId
uint32_t setDetUnitId(const G4Step *step) override
Definition: AHCalSD.cc:87
AHCalParameters::kRow_
static constexpr int kRow_
Definition: AHCalParameters.h:28
CaloSD::getAttenuation
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:517
MeV
const double MeV
AHCalSD::eminHit
double eminHit
Definition: AHCalSD.cc:43
AHCalSD::getEnergyDeposit
double getEnergyDeposit(const G4Step *) override
Definition: AHCalSD.cc:71
AHCalSD::betaThr
double betaThr
Definition: AHCalSD.cc:42
AHCalSD::useBirk
bool useBirk
Definition: AHCalSD.cc:41
heavyIonCSV_trainingSettings.idx
idx
Definition: heavyIonCSV_trainingSettings.py:5
DetId
Definition: DetId.h:17
MakerMacros.h
AHCalDetId::depth
int depth() const
get the layer number
Definition: AHCalDetId.cc:43
CaloG4Hit::getEnergyDeposit
double getEnergyDeposit() const
Definition: CaloG4Hit.h:77
AHCalSD::birk1
double birk1
Definition: AHCalSD.cc:42
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
AHCalDetId::icol
int icol() const
get the column number
Definition: AHCalDetId.cc:36
gen
Definition: PythiaDecays.h:13
LEDCalibrationChannels.depth
depth
Definition: LEDCalibrationChannels.py:65
AHCalSD::birk3
double birk3
Definition: AHCalSD.cc:42
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
AHCalDetId
Definition: AHCalDetId.h:13
DEFINE_SENSITIVEDETECTOR
#define DEFINE_SENSITIVEDETECTOR(type)
Definition: SensitiveDetectorPluginFactory.h:13
AHCalParameters::kSignRowColumn_
static constexpr int kSignRowColumn_
Definition: AHCalParameters.h:31
edm::ParameterSet
Definition: ParameterSet.h:47
CaloSD.h
ParameterSet
Definition: Functions.h:16
HcalDetId.h
ModuleDef.h
edm::EventSetup
Definition: EventSetup.h:57
TrackInformation.h
CaloG4Hit
Definition: CaloG4Hit.h:32
AHCalParameters::kColumn_
static constexpr int kColumn_
Constants used.
Definition: AHCalParameters.h:27
AHcalSensitiveDetector
AHCalSD AHcalSensitiveDetector
Definition: AHCalSD.cc:134
HcalSubdetector
HcalSubdetector
Definition: HcalAssistant.h:31
AHCalSD::filterHit
bool filterHit(CaloG4Hit *, double) override
Definition: AHCalSD.cc:126
DetId::rawId
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
AHCalSD::AHCalSD
AHCalSD(const std::string &, const edm::EventSetup &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: AHCalSD.cc:46
AHCalParameters.h
DetId.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
relval_steps.gen
def gen(fragment, howMuch)
Production test section ####.
Definition: relval_steps.py:509
AlignmentPI::index
index
Definition: AlignmentPayloadInspectorHelper.h:46
AHCalSD
Definition: AHCalSD.cc:25
SensitiveDetectorPluginFactory.h
ntuplemaker.time
time
Definition: ntuplemaker.py:310
AHCalParameters::kRowColumn_
static constexpr int kRowColumn_
Definition: AHCalParameters.h:30
TauDecayModes.dec
dec
Definition: TauDecayModes.py:143
AHCalDetId.h
weight
Definition: weight.py:1
CaloSD
Definition: CaloSD.h:38
g
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
RemoveAddSevLevel.flag
flag
Definition: RemoveAddSevLevel.py:116
AHCalParameters::kSign_
static constexpr int kSign_
Definition: AHCalParameters.h:29