CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
EcalTBH4BeamSD.cc
Go to the documentation of this file.
1 // File: EcalTBH4BeamSD.cc
3 // Description: Sensitive Detector class for electromagnetic calorimeters
13 
15 
16 #include "G4Step.hh"
17 #include "G4Track.hh"
18 #include "G4VProcess.hh"
19 
20 #include "G4SystemOfUnits.hh"
21 
23  const SensitiveDetectorCatalog &clg,
24  edm::ParameterSet const &p,
25  const SimTrackManager *manager)
26  : CaloSD(name, clg, p, manager), numberingScheme(nullptr) {
27  edm::ParameterSet m_EcalTBH4BeamSD = p.getParameter<edm::ParameterSet>("EcalTBH4BeamSD");
28  useBirk = m_EcalTBH4BeamSD.getParameter<bool>("UseBirkLaw");
29  birk1 = m_EcalTBH4BeamSD.getParameter<double>("BirkC1") * (g / (MeV * cm2));
30  birk2 = m_EcalTBH4BeamSD.getParameter<double>("BirkC2");
31  birk3 = m_EcalTBH4BeamSD.getParameter<double>("BirkC3");
32 
33  EcalNumberingScheme *scheme = nullptr;
34  if (name == "EcalTBH4BeamHits") {
35  scheme = dynamic_cast<EcalNumberingScheme *>(new EcalHodoscopeNumberingScheme());
36  } else {
37  edm::LogWarning("EcalTBSim") << "EcalTBH4BeamSD: ReadoutName not supported\n";
38  }
39 
40  if (scheme)
41  setNumberingScheme(scheme);
42  edm::LogInfo("EcalTBSim") << "Constructing a EcalTBH4BeamSD with name " << GetName();
43  edm::LogInfo("EcalTBSim") << "EcalTBH4BeamSD:: Use of Birks law is set to " << useBirk
44  << " with three constants kB = " << birk1 << ", C1 = " << birk2
45  << ", C2 = " << birk3;
46 }
47 
49  if (numberingScheme)
50  delete numberingScheme;
51 }
52 
53 double EcalTBH4BeamSD::getEnergyDeposit(const G4Step *aStep) {
54  // take into account light collection curve for crystals
55  double weight = 1.;
56  if (useBirk)
57  weight *= getAttenuation(aStep, birk1, birk2, birk3);
58  double edep = aStep->GetTotalEnergyDeposit() * weight;
59  LogDebug("EcalTBSim") << "EcalTBH4BeamSD:: " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()
60  << " Light Collection Efficiency " << weight << " Weighted Energy Deposit " << edep / MeV
61  << " MeV";
62  return edep;
63 }
64 
65 uint32_t EcalTBH4BeamSD::setDetUnitId(const G4Step *aStep) {
66  getBaseNumber(aStep);
67  return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(theBaseNumber));
68 }
69 
71  if (scheme != nullptr) {
72  edm::LogInfo("EcalTBSim") << "EcalTBH4BeamSD: updates numbering scheme for " << GetName() << "\n";
73  if (numberingScheme)
74  delete numberingScheme;
75  numberingScheme = scheme;
76  }
77 }
78 
79 void EcalTBH4BeamSD::getBaseNumber(const G4Step *aStep) {
81  const G4VTouchable *touch = aStep->GetPreStepPoint()->GetTouchable();
82  int theSize = touch->GetHistoryDepth() + 1;
83  if (theBaseNumber.getCapacity() < theSize)
84  theBaseNumber.setSize(theSize);
85  // Get name and copy numbers
86  if (theSize > 1) {
87  for (int ii = 0; ii < theSize; ii++) {
88  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
89  LogDebug("EcalTBSim") << "EcalTBH4BeamSD::getBaseNumber(): Adding level " << ii << ": "
90  << touch->GetVolume(ii)->GetName() << "[" << touch->GetReplicaNumber(ii) << "]";
91  }
92  }
93 }
Definition: CaloSD.h:40
int ii
Definition: cuy.py:589
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
void addLevel(const std::string &name, const int &copyNumber)
double getEnergyDeposit(const G4Step *) override
uint32_t setDetUnitId(const G4Step *step) override
void setNumberingScheme(EcalNumberingScheme *scheme)
EcalTBH4BeamSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
EcalNumberingScheme * numberingScheme
Log< level::Info, false > LogInfo
EcalBaseNumber theBaseNumber
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
~EcalTBH4BeamSD() override
void getBaseNumber(const G4Step *aStep)
static TrackerG4SimHitNumberingScheme & numberingScheme(const GeometricDet &det)
Log< level::Warning, false > LogWarning
int weight
Definition: histoStyle.py:51
void setSize(const int &size)
virtual uint32_t getUnitID(const EcalBaseNumber &baseNumber) const =0
#define LogDebug(id)