Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:21:18

0001 ///////////////////////////////////////////////////////////////////////////////
0002 // File: EcalTBH4BeamSD.cc
0003 // Description: Sensitive Detector class for electromagnetic calorimeters
0004 ///////////////////////////////////////////////////////////////////////////////
0005 #include "DetectorDescription/Core/interface/DDFilter.h"
0006 #include "DetectorDescription/Core/interface/DDFilteredView.h"
0007 #include "DetectorDescription/Core/interface/DDSolid.h"
0008 #include "DetectorDescription/Core/interface/DDSplit.h"
0009 #include "DetectorDescription/Core/interface/DDValue.h"
0010 #include "Geometry/EcalCommonData/interface/EcalBaseNumber.h"
0011 #include "Geometry/EcalTestBeam/interface/EcalHodoscopeNumberingScheme.h"
0012 #include "SimG4CMS/EcalTestBeam/interface/EcalTBH4BeamSD.h"
0013 
0014 #include "Geometry/EcalCommonData/interface/EcalBaseNumber.h"
0015 
0016 #include "G4Step.hh"
0017 #include "G4Track.hh"
0018 #include "G4VProcess.hh"
0019 
0020 #include <CLHEP/Units/SystemOfUnits.h>
0021 
0022 EcalTBH4BeamSD::EcalTBH4BeamSD(const std::string &name,
0023                                const SensitiveDetectorCatalog &clg,
0024                                edm::ParameterSet const &p,
0025                                const SimTrackManager *manager)
0026     : CaloSD(name, clg, p, manager), numberingScheme(nullptr) {
0027   edm::ParameterSet m_EcalTBH4BeamSD = p.getParameter<edm::ParameterSet>("EcalTBH4BeamSD");
0028   useBirk = m_EcalTBH4BeamSD.getParameter<bool>("UseBirkLaw");
0029   birk1 = m_EcalTBH4BeamSD.getParameter<double>("BirkC1") * (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
0030   birk2 = m_EcalTBH4BeamSD.getParameter<double>("BirkC2");
0031   birk3 = m_EcalTBH4BeamSD.getParameter<double>("BirkC3");
0032 
0033   EcalNumberingScheme *scheme = nullptr;
0034   if (name == "EcalTBH4BeamHits") {
0035     scheme = dynamic_cast<EcalNumberingScheme *>(new EcalHodoscopeNumberingScheme());
0036   } else {
0037     edm::LogWarning("EcalTBSim") << "EcalTBH4BeamSD: ReadoutName not supported\n";
0038   }
0039 
0040   if (scheme)
0041     setNumberingScheme(scheme);
0042   edm::LogInfo("EcalTBSim") << "Constructing a EcalTBH4BeamSD  with name " << GetName();
0043   edm::LogInfo("EcalTBSim") << "EcalTBH4BeamSD:: Use of Birks law is set to  " << useBirk
0044                             << "        with three constants kB = " << birk1 << ", C1 = " << birk2
0045                             << ", C2 = " << birk3;
0046 }
0047 
0048 EcalTBH4BeamSD::~EcalTBH4BeamSD() {
0049   if (numberingScheme)
0050     delete numberingScheme;
0051 }
0052 
0053 double EcalTBH4BeamSD::getEnergyDeposit(const G4Step *aStep) {
0054   // take into account light collection curve for crystals
0055   double weight = 1.;
0056   if (useBirk)
0057     weight *= getAttenuation(aStep, birk1, birk2, birk3);
0058   double edep = aStep->GetTotalEnergyDeposit() * weight;
0059   LogDebug("EcalTBSim") << "EcalTBH4BeamSD:: " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetName()
0060                         << " Light Collection Efficiency " << weight << " Weighted Energy Deposit " << edep / CLHEP::MeV
0061                         << " MeV";
0062   return edep;
0063 }
0064 
0065 uint32_t EcalTBH4BeamSD::setDetUnitId(const G4Step *aStep) {
0066   getBaseNumber(aStep);
0067   return (numberingScheme == nullptr ? 0 : numberingScheme->getUnitID(theBaseNumber));
0068 }
0069 
0070 void EcalTBH4BeamSD::setNumberingScheme(EcalNumberingScheme *scheme) {
0071   if (scheme != nullptr) {
0072     edm::LogInfo("EcalTBSim") << "EcalTBH4BeamSD: updates numbering scheme for " << GetName() << "\n";
0073     if (numberingScheme)
0074       delete numberingScheme;
0075     numberingScheme = scheme;
0076   }
0077 }
0078 
0079 void EcalTBH4BeamSD::getBaseNumber(const G4Step *aStep) {
0080   theBaseNumber.reset();
0081   const G4VTouchable *touch = aStep->GetPreStepPoint()->GetTouchable();
0082   int theSize = touch->GetHistoryDepth() + 1;
0083   if (theBaseNumber.getCapacity() < theSize)
0084     theBaseNumber.setSize(theSize);
0085   // Get name and copy numbers
0086   if (theSize > 1) {
0087     for (int ii = 0; ii < theSize; ii++) {
0088       theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(), touch->GetReplicaNumber(ii));
0089       LogDebug("EcalTBSim") << "EcalTBH4BeamSD::getBaseNumber(): Adding level " << ii << ": "
0090                             << touch->GetVolume(ii)->GetName() << "[" << touch->GetReplicaNumber(ii) << "]";
0091     }
0092   }
0093 }