CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TotemT2ScintSD.cc
Go to the documentation of this file.
5 
6 #include "G4SDManager.hh"
7 #include "G4Step.hh"
8 #include "G4Track.hh"
9 #include "G4VProcess.hh"
10 #include "G4ios.hh"
11 #include "G4Cerenkov.hh"
12 #include "G4ParticleTable.hh"
13 #include "CLHEP/Units/GlobalSystemOfUnits.h"
14 #include "CLHEP/Units/GlobalPhysicalConstants.h"
15 #include "Randomize.hh"
16 #include "G4Poisson.hh"
17 
18 //#define EDM_ML_DEBUG
19 
21  const SensitiveDetectorCatalog& clg,
22  edm::ParameterSet const& p,
23  const SimTrackManager* manager)
24  : CaloSD(name,
25  clg,
26  p,
27  manager,
28  (float)(p.getParameter<edm::ParameterSet>("TotemT2ScintSD").getParameter<double>("TimeSliceUnit")),
29  p.getParameter<edm::ParameterSet>("TotemT2ScintSD").getParameter<bool>("IgnoreTrackID")) {
30  edm::ParameterSet m_T2SD = p.getParameter<edm::ParameterSet>("TotemT2ScintSD");
31  useBirk_ = m_T2SD.getParameter<bool>("UseBirkLaw");
32  birk1_ = m_T2SD.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
33  birk2_ = m_T2SD.getParameter<double>("BirkC2");
34  birk3_ = m_T2SD.getParameter<double>("BirkC3");
36 
37  edm::LogVerbatim("ForwardSim") << "***************************************************\n"
38  << "* *\n"
39  << "* Constructing a TotemT2ScintSD with name " << name << " *\n"
40  << "* *\n"
41  << "***************************************************";
42 
43  edm::LogVerbatim("ForwardSim") << "\nUse of Birks law is set to " << useBirk_
44  << " with three constants kB = " << birk1_ << ", C1 = " << birk2_
45  << ", C2 = " << birk3_;
46 }
47 
48 uint32_t TotemT2ScintSD::setDetUnitId(const G4Step* aStep) {
49  auto const prePoint = aStep->GetPreStepPoint();
50  auto const touch = prePoint->GetTouchable();
51 
52  int iphi = (touch->GetReplicaNumber(0)) % 10;
53  int lay = (touch->GetReplicaNumber(0) / 10) % 10 + 1;
54  int zside = (((touch->GetReplicaNumber(1)) == 1) ? 1 : -1);
55 
56  return setDetUnitId(zside, lay, iphi);
57 }
58 
60  if (scheme != nullptr) {
61  edm::LogVerbatim("ForwardSim") << "TotemT2ScintSD: updates numbering scheme for " << GetName();
62  numberingScheme.reset(scheme);
63  }
64 }
65 
66 double TotemT2ScintSD::getEnergyDeposit(const G4Step* aStep) {
67  double destep = aStep->GetTotalEnergyDeposit();
68  double weight = ((useBirk_) ? getAttenuation(aStep, birk1_, birk2_, birk3_) : 1.0);
69  double edep = weight * destep;
70 #ifdef EDM_ML_DEBUG
71  edm::LogVerbatim("ForwardSim") << "TotemT2ScintSD: edep= " << destep << ":" << weight << ":" << edep;
72 #endif
73  return edep;
74 }
75 
76 uint32_t TotemT2ScintSD::setDetUnitId(const int& zside, const int& lay, const int& iphi) {
77  uint32_t id = ((numberingScheme.get()) ? numberingScheme->packID(zside, lay, iphi) : 0);
78 #ifdef EDM_ML_DEBUG
79  edm::LogVerbatim("ForwardSim") << "TotemT2ScintSD: zside " << zside << " layer " << lay << " phi " << iphi << " ID "
80  << std::hex << id << std::dec;
81 #endif
82  return id;
83 }
Log< level::Info, true > LogVerbatim
uint32_t setDetUnitId(const G4Step *step) override
double getEnergyDeposit(const G4Step *) override
uint16_t *__restrict__ id
Definition: CaloSD.h:40
int zside(DetId const &)
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
const double MeV
TotemT2ScintSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
void setNumberingScheme(TotemT2ScintNumberingScheme *scheme)
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:657
std::unique_ptr< TotemT2ScintNumberingScheme > numberingScheme
int weight
Definition: histoStyle.py:51