CMS 3D CMS Logo

AHCalSD.cc
Go to the documentation of this file.
6 
7 #include "G4LogicalVolumeStore.hh"
8 #include "G4LogicalVolume.hh"
9 #include "G4Step.hh"
10 #include "G4Track.hh"
11 #include "G4ParticleTable.hh"
12 #include "G4VProcess.hh"
13 
14 #include "G4SystemOfUnits.hh"
15 #include "G4PhysicalConstants.hh"
16 
17 #include <iomanip>
18 //#define EDM_ML_DEBUG
19 
21  const SensitiveDetectorCatalog & clg,
22  edm::ParameterSet const & p, const SimTrackManager* manager) :
23  CaloSD(name, cpv, clg, p, manager,
24  (float)(p.getParameter<edm::ParameterSet>("AHCalSD").getParameter<double>("TimeSliceUnit")),
25  p.getParameter<edm::ParameterSet>("AHCalSD").getParameter<bool>("IgnoreTrackID")) {
26 
28  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
29  birk1 = m_HC.getParameter<double>("BirkC1")*(CLHEP::g/(CLHEP::MeV*CLHEP::cm2));
30  birk2 = m_HC.getParameter<double>("BirkC2");
31  birk3 = m_HC.getParameter<double>("BirkC3");
32  eminHit = m_HC.getParameter<double>("EminHit")*CLHEP::MeV;
33 
34  edm::LogInfo("HcalSim") << "AHCalSD:: Use of Birks law is set to "
35  << useBirk << " with three constants kB = "
36  << birk1 << ", C1 = " << birk2 << ", C2 = " << birk3
37  << "\nAHCalSD:: Threshold for storing hits: "
38  << eminHit << std::endl;
39 }
40 
41 double AHCalSD::getEnergyDeposit(const G4Step* aStep) {
42 
43  double destep = aStep->GetTotalEnergyDeposit();
44  double wt2 = aStep->GetTrack()->GetWeight();
45  double weight = (wt2 > 0.0) ? wt2 : 1.0;
46 #ifdef EDM_ML_DEBUG
47  double weight0 = weight;
48 #endif
49  if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
50 #ifdef EDM_ML_DEBUG
51  edm::LogInfo("HcalSim") << "AHCalSD: weight " << weight0 << " " << weight
52  << std::endl;
53 #endif
54  double edep = weight*destep;
55  return edep;
56 }
57 
58 uint32_t AHCalSD::setDetUnitId(const G4Step * aStep) {
59 
60  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
61  const G4VTouchable* touch = preStepPoint->GetTouchable();
62 
63  int depth = (touch->GetReplicaNumber(1));
64  int incol = ((touch->GetReplicaNumber(0))%10);
65  int inrow = ((touch->GetReplicaNumber(0))/10)%10;
66  int jncol = ((touch->GetReplicaNumber(0))/100)%10;
67  int jnrow = ((touch->GetReplicaNumber(0))/1000)%10;
68  int col = (jncol == 0) ? incol : -incol;
69  int row = (jnrow == 0) ? inrow : -inrow;
70  uint32_t index = AHCalDetId(row,col,depth).rawId();
71 #ifdef EDM_ML_DEBUG
72  edm::LogInfo("HcalSim") << "AHCalSD: det = " << HcalOther
73  << " depth = " << depth << " row = " << row
74  << " column = " << col << " packed index = 0x"
75  << std::hex << index << std::dec << std::endl;
76  bool flag = unpackIndex(index, row, col, depth);
77  edm::LogInfo("HcalSim") << "Results from unpacker for 0x" << std::hex
78  << index << std::dec << " Flag " << flag << " Row "
79  << row << " Col " << col << " Depth " << depth
80  << std::endl;
81 #endif
82  return index;
83 }
84 
85 bool AHCalSD::unpackIndex(const uint32_t& idx, int& row, int& col, int& depth) {
86 
87  DetId gen(idx);
88  HcalSubdetector subdet = (HcalSubdetector(gen.subdetId()));
89  bool rcode = (gen.det()==DetId::Hcal && subdet!=HcalOther);
90  row = col = depth = 0;
91  if (rcode) {
92  row = AHCalDetId(idx).irow();
93  col = AHCalDetId(idx).icol();
94  depth = AHCalDetId(idx).depth();
95  }
96 #ifdef EDM_ML_DEBUG
97  edm::LogInfo("HcalSim") << "AHCalSD: packed index = 0x" << std::hex << idx
98  << std::dec << " Row " << row << " Column " << col
99  << " Depth " << depth << " OK " << rcode << std::endl;
100 #endif
101  return rcode;
102 }
103 
104 bool AHCalSD::filterHit(CaloG4Hit* aHit, double time) {
105  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit));
106 }
T getParameter(std::string const &) const
double birk1
Definition: AHCalSD.h:29
Definition: CaloSD.h:37
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:47
double birk3
Definition: AHCalSD.h:29
int irow() const
get the row number
Definition: AHCalDetId.cc:33
Definition: weight.py:1
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
bool filterHit(CaloG4Hit *, double) override
Definition: AHCalSD.cc:104
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
int depth() const
Definition: AHCalDetId.cc:45
const double MeV
int icol() const
get the column number
Definition: AHCalDetId.cc:39
constexpr int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:41
HcalSubdetector
Definition: HcalAssistant.h:31
double birk2
Definition: AHCalSD.h:29
double getEnergyDeposit(const G4Step *) override
Definition: AHCalSD.cc:41
uint32_t setDetUnitId(const G4Step *step) override
Definition: AHCalSD.cc:58
def gen(fragment, howMuch)
Production test section ####.
AHCalSD(const std::string &, const DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
Definition: AHCalSD.cc:20
double tmaxHit
Definition: CaloSD.h:132
Definition: DetId.h:18
double eminHit
Definition: AHCalSD.h:30
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:462
HLT enums.
col
Definition: cuy.py:1009
bool useBirk
Definition: AHCalSD.h:28
bool unpackIndex(const uint32_t &idx, int &row, int &col, int &depth)
Definition: AHCalSD.cc:85
double getEnergyDeposit() const
Definition: CaloG4Hit.h:81
constexpr Detector det() const
get the detector field from this detid
Definition: DetId.h:39