CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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 
24  edm::ParameterSet const & p,
25  const SimTrackManager* manager) :
26  CaloSD(name, cpv, clg, p, manager), numberingScheme(0) {
27 
28  edm::ParameterSet m_EcalTBH4BeamSD = p.getParameter<edm::ParameterSet>("EcalTBH4BeamSD");
29  useBirk= m_EcalTBH4BeamSD.getParameter<bool>("UseBirkLaw");
30  birk1 = m_EcalTBH4BeamSD.getParameter<double>("BirkC1")*(g/(MeV*cm2));
31  birk2 = m_EcalTBH4BeamSD.getParameter<double>("BirkC2");
32  birk3 = m_EcalTBH4BeamSD.getParameter<double>("BirkC3");
33 
34  EcalNumberingScheme* scheme=0;
35  if (name == "EcalTBH4BeamHits") {
36  scheme = dynamic_cast<EcalNumberingScheme*>(new EcalHodoscopeNumberingScheme());
37  }
38  else {edm::LogWarning("EcalTBSim") << "EcalTBH4BeamSD: ReadoutName not supported\n";}
39 
40  if (scheme) setNumberingScheme(scheme);
41  edm::LogInfo("EcalTBSim") << "Constructing a EcalTBH4BeamSD with name "
42  << GetName();
43  edm::LogInfo("EcalTBSim") << "EcalTBH4BeamSD:: Use of Birks law is set to "
44  << useBirk << " with three constants kB = "
45  << birk1 << ", C1 = " << birk2 << ", C2 = "
46  << birk3;
47 }
48 
50  if (numberingScheme) delete numberingScheme;
51 }
52 
53 double EcalTBH4BeamSD::getEnergyDeposit(G4Step * aStep) {
54 
55  if (aStep == NULL) {
56  return 0;
57  } else {
58  preStepPoint = aStep->GetPreStepPoint();
59  G4String nameVolume = preStepPoint->GetPhysicalVolume()->GetName();
60 
61  // take into account light collection curve for crystals
62  double weight = 1.;
63  if (useBirk) weight *= getAttenuation(aStep, birk1, birk2, birk3);
64  double edep = aStep->GetTotalEnergyDeposit() * weight;
65  LogDebug("EcalTBSim") << "EcalTBH4BeamSD:: " << nameVolume
66  <<" Light Collection Efficiency " << weight
67  << " Weighted Energy Deposit " << edep/MeV << " MeV";
68  return edep;
69  }
70 }
71 
72 uint32_t EcalTBH4BeamSD::setDetUnitId(G4Step * aStep) {
73  getBaseNumber(aStep);
75 }
76 
78  if (scheme != 0) {
79  edm::LogInfo("EcalTBSim") << "EcalTBH4BeamSD: updates numbering scheme for "
80  << GetName() << "\n";
81  if (numberingScheme) delete numberingScheme;
82  numberingScheme = scheme;
83  }
84 }
85 
86 
87 void EcalTBH4BeamSD::getBaseNumber(const G4Step* aStep) {
88 
90  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
91  int theSize = touch->GetHistoryDepth()+1;
92  if ( theBaseNumber.getCapacity() < theSize ) theBaseNumber.setSize(theSize);
93  //Get name and copy numbers
94  if ( theSize > 1 ) {
95  for (int ii = 0; ii < theSize ; ii++) {
96  theBaseNumber.addLevel(touch->GetVolume(ii)->GetName(),touch->GetReplicaNumber(ii));
97  LogDebug("EcalTBSim") << "EcalTBH4BeamSD::getBaseNumber(): Adding level " << ii
98  << ": " << touch->GetVolume(ii)->GetName() << "["
99  << touch->GetReplicaNumber(ii) << "]";
100  }
101  }
102 }
#define LogDebug(id)
T getParameter(std::string const &) const
Definition: CaloSD.h:42
virtual uint32_t setDetUnitId(G4Step *step)
#define NULL
Definition: scimark2.h:8
int ii
Definition: cuy.py:588
type of data representation of DDCompactView
Definition: DDCompactView.h:77
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
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
const double MeV
void addLevel(const std::string &name, const int &copyNumber)
virtual ~EcalTBH4BeamSD()
void setNumberingScheme(EcalNumberingScheme *scheme)
double getAttenuation(G4Step *aStep, double birk1, double birk2, double birk3)
Definition: CaloSD.cc:465
EcalNumberingScheme * numberingScheme
G4StepPoint * preStepPoint
Definition: CaloSD.h:120
EcalBaseNumber theBaseNumber
EcalTBH4BeamSD(G4String, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
void getBaseNumber(const G4Step *aStep)
virtual double getEnergyDeposit(G4Step *)
int weight
Definition: histoStyle.py:50
void setSize(const int &size)
virtual uint32_t getUnitID(const EcalBaseNumber &baseNumber) const =0