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 
42 
43 double AHCalSD::getEnergyDeposit(G4Step* aStep) {
44 
45  double destep = aStep->GetTotalEnergyDeposit();
46  double wt2 = aStep->GetTrack()->GetWeight();
47  double weight = (wt2 > 0.0) ? wt2 : 1.0;
48 #ifdef EDM_ML_DEBUG
49  double weight0 = weight;
50 #endif
51  if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
52 #ifdef EDM_ML_DEBUG
53  edm::LogInfo("HcalSim") << "AHCalSD: weight " << weight0 << " " << weight
54  << std::endl;
55 #endif
56  double edep = weight*destep;
57  return edep;
58 }
59 
60 uint32_t AHCalSD::setDetUnitId(const G4Step * aStep) {
61 
62  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
63  const G4VTouchable* touch = preStepPoint->GetTouchable();
64 
65  int depth = (touch->GetReplicaNumber(1));
66  int incol = ((touch->GetReplicaNumber(0))%10);
67  int inrow = ((touch->GetReplicaNumber(0))/10)%10;
68  int jncol = ((touch->GetReplicaNumber(0))/100)%10;
69  int jnrow = ((touch->GetReplicaNumber(0))/1000)%10;
70  int col = (jncol == 0) ? incol : -incol;
71  int row = (jnrow == 0) ? inrow : -inrow;
72  uint32_t index = AHCalDetId(row,col,depth).rawId();
73 #ifdef EDM_ML_DEBUG
74  edm::LogInfo("HcalSim") << "AHCalSD: det = " << HcalOther
75  << " depth = " << depth << " row = " << row
76  << " column = " << col << " packed index = 0x"
77  << std::hex << index << std::dec << std::endl;
78  bool flag = unpackIndex(index, row, col, depth);
79  edm::LogInfo("HcalSim") << "Results from unpacker for 0x" << std::hex
80  << index << std::dec << " Flag " << flag << " Row "
81  << row << " Col " << col << " Depth " << depth
82  << std::endl;
83 #endif
84  return index;
85 }
86 
87 bool AHCalSD::unpackIndex(const uint32_t& idx, int& row, int& col, int& depth) {
88 
89  DetId gen(idx);
90  HcalSubdetector subdet = (HcalSubdetector(gen.subdetId()));
91  bool rcode = (gen.det()==DetId::Hcal && subdet!=HcalOther);
92  row = col = depth = 0;
93  if (rcode) {
94  row = AHCalDetId(idx).irow();
95  col = AHCalDetId(idx).icol();
96  depth = AHCalDetId(idx).depth();
97  }
98 #ifdef EDM_ML_DEBUG
99  edm::LogInfo("HcalSim") << "AHCalSD: packed index = 0x" << std::hex << idx
100  << std::dec << " Row " << row << " Column " << col
101  << " Depth " << depth << " OK " << rcode << std::endl;
102 #endif
103  return rcode;
104 }
105 
106 bool AHCalSD::filterHit(CaloG4Hit* aHit, double time) {
107  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > eminHit));
108 }
T getParameter(std::string const &) const
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:439
double birk1
Definition: AHCalSD.h:30
Definition: CaloSD.h:42
double birk3
Definition: AHCalSD.h:30
int irow() const
get the row number
Definition: AHCalDetId.cc:33
Definition: weight.py:1
type of data representation of DDCompactView
Definition: DDCompactView.h:90
bool filterHit(CaloG4Hit *, double) override
Definition: AHCalSD.cc:106
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
uint32_t rawId() const
get the raw id
Definition: DetId.h:43
const double MeV
int icol() const
get the column number
Definition: AHCalDetId.cc:39
HcalSubdetector
Definition: HcalAssistant.h:31
double birk2
Definition: AHCalSD.h:30
~AHCalSD() override
Definition: AHCalSD.cc:41
uint32_t setDetUnitId(const G4Step *step) override
Definition: AHCalSD.cc:60
int subdetId() const
get the contents of the subdetector field (not cast into any detector&#39;s numbering enum) ...
Definition: DetId.h:37
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 getEnergyDeposit(G4Step *) override
Definition: AHCalSD.cc:43
double tmaxHit
Definition: CaloSD.h:124
G4StepPoint * preStepPoint
Definition: CaloSD.h:121
Definition: DetId.h:18
double eminHit
Definition: AHCalSD.h:31
HLT enums.
col
Definition: cuy.py:1008
bool useBirk
Definition: AHCalSD.h:29
Detector det() const
get the detector field from this detid
Definition: DetId.h:35
bool unpackIndex(const uint32_t &idx, int &row, int &col, int &depth)
Definition: AHCalSD.cc:87
double getEnergyDeposit() const
Definition: CaloG4Hit.h:81